]> git.donarmstrong.com Git - mothur.git/blob - shhhercommand.cpp
mods to shhh.seqs command
[mothur.git] / shhhercommand.cpp
1 /*
2  *  shhher.cpp
3  *  Mothur
4  *
5  *  Created by Pat Schloss on 12/27/10.
6  *  Copyright 2010 Schloss Lab. All rights reserved.
7  *
8  */
9
10 #include "shhhercommand.h"
11
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"
19 #include <cfloat>
20
21 //**********************************************************************************************************************
22
23 #define NUMBINS 1000
24 #define HOMOPS 10
25 #define MIN_COUNT 0.1
26 #define MIN_WEIGHT 0.1
27 #define MIN_TAU 0.0001
28 #define MIN_ITER 10
29
30 //**********************************************************************************************************************
31
32 vector<string> ShhherCommand::getValidParameters(){     
33         try {
34                 string Array[] =  {     
35                         "file", "flow", "lookup", "cutoff", "sigma", "outputdir","inputdir", "processors"       
36                 };
37                 
38                 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
39                 return myArray;
40         }
41         catch(exception& e) {
42                 m->errorOut(e, "ShhherCommand", "getValidParameters");
43                 exit(1);
44         }
45 }
46
47 //**********************************************************************************************************************
48
49 ShhherCommand::ShhherCommand(){ 
50         try {
51                 abort = true; calledHelp = true;
52                 
53                 //initialize outputTypes
54                 vector<string> tempOutNames;
55                 outputTypes["pn.dist"] = tempOutNames;
56
57         }
58         catch(exception& e) {
59                 m->errorOut(e, "ShhherCommand", "ShhherCommand");
60                 exit(1);
61         }
62 }
63
64 //**********************************************************************************************************************
65
66 vector<string> ShhherCommand::getRequiredParameters(){  
67         try {
68                 string Array[] =  {"flow"};
69                 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
70                 return myArray;
71         }
72         catch(exception& e) {
73                 m->errorOut(e, "ShhherCommand", "getRequiredParameters");
74                 exit(1);
75         }
76 }
77
78 //**********************************************************************************************************************
79
80 vector<string> ShhherCommand::getRequiredFiles(){       
81         try {
82                 vector<string> myArray;
83                 return myArray;
84         }
85         catch(exception& e) {
86                 m->errorOut(e, "ShhherCommand", "getRequiredFiles");
87                 exit(1);
88         }
89 }
90
91 //**********************************************************************************************************************
92
93 ShhherCommand::ShhherCommand(string option) {
94         try {
95
96 #ifdef USE_MPI
97                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
98                 MPI_Comm_size(MPI_COMM_WORLD, &ncpus);
99
100                 if(pid == 0){
101 #endif
102                 
103                 
104                 abort = false; calledHelp = false;   
105                 
106                 
107                 //allow user to run help
108                 if(option == "help") { help(); abort = true; calledHelp = true; }
109                 
110                 else {
111                         
112                         //valid paramters for this command
113                         string AlignArray[] =  {
114                                 "file", "flow", "lookup", "cutoff", "sigma", "outputdir","inputdir", "processors"       
115                         };
116                         
117                         vector<string> myArray (AlignArray, AlignArray+(sizeof(AlignArray)/sizeof(string)));
118                         
119                         OptionParser parser(option);
120                         map<string,string> parameters = parser.getParameters();
121                         
122                         ValidParameters validParameter;
123                         map<string,string>::iterator it;
124                         
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;  }
128                         }
129                         
130                         //initialize outputTypes
131                         vector<string> tempOutNames;
132                         outputTypes["pn.dist"] = tempOutNames;
133                         //                      outputTypes["fasta"] = tempOutNames;
134                         
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 = "";          }
138                         else {
139                                 string path;
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;             }
146                                 }
147                                 
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;           }
154                                 }
155
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;             }
162                                 }
163                         }
164                         
165                         
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();
172                                 abort = true; 
173                         }
174                         else if (flowFileName == "not open" || flowFilesFileName == "not open") { abort = true; }
175                         
176                         if(flowFileName != "not found"){        compositeFASTAFileName = "";    }
177                         else{
178                                 compositeFASTAFileName = flowFilesFileName.substr(0, flowFilesFileName.length()-10) + "pn.fasta";
179                                 ofstream temp;
180                                 m->openOutputFile(compositeFASTAFileName, temp);
181                                 temp.close();
182                         }
183                         
184                         //if the user changes the output directory command factory will send this info to us in the output parameter 
185                         outputDir = validParameter.validFile(parameters, "outputdir", false);           if (outputDir == "not found"){  
186                                 outputDir = ""; 
187                                 outputDir += m->hasPath(flowFileName); //if user entered a file with a path then preserve it    
188                         }
189                         
190                         
191                         //check for optional parameter and set defaults
192                         // ...at some point should added some additional type checking...
193                         string temp;
194                         temp = validParameter.validFile(parameters, "lookup", true);
195                         if (temp == "not found")        {       lookupFileName = "LookUp_Titanium.pat"; }
196                         else if(temp == "not open")     {       abort = true;                   } 
197                         else                                            {       lookupFileName = temp;  }
198                         
199                         temp = validParameter.validFile(parameters, "processors", false);if (temp == "not found"){      temp = "1";                     }
200                         convert(temp, processors); 
201
202                         temp = validParameter.validFile(parameters, "cutoff", false);   if (temp == "not found"){       temp = "0.01";          }
203                         convert(temp, cutoff); 
204                         
205                         temp = validParameter.validFile(parameters, "mindelta", false); if (temp == "not found"){       temp = "0.000001";      }
206                         convert(temp, minDelta); 
207
208                         temp = validParameter.validFile(parameters, "maxiter", false);  if (temp == "not found"){       temp = "1000";          }
209                         convert(temp, maxIters); 
210
211                         temp = validParameter.validFile(parameters, "sigma", false);if (temp == "not found")    {       temp = "60";            }
212                         convert(temp, sigma); 
213                         
214                         globaldata = GlobalData::getInstance();
215                 }
216                         
217 #ifdef USE_MPI
218                 }                               
219 #endif
220                                 
221         }
222         catch(exception& e) {
223                 m->errorOut(e, "ShhherCommand", "ShhherCommand");
224                 exit(1);
225         }
226 }
227
228 //**********************************************************************************************************************
229
230 ShhherCommand::~ShhherCommand(){}
231
232 //**********************************************************************************************************************
233
234 void ShhherCommand::help(){
235         try {
236                 m->mothurOut("The shhher command reads a file containing flowgrams and creates a file of corrected sequences.\n");
237         }
238         catch(exception& e) {
239                 m->errorOut(e, "ShhherCommand", "help");
240                 exit(1);
241         }
242 }
243
244 //**********************************************************************************************************************
245 #ifdef USE_MPI
246 int ShhherCommand::execute(){
247         try {
248                 if (abort == true) { if (calledHelp) { return 0; }  return 2;   }
249                 
250                 int tag = 1976;
251                 MPI_Status status; 
252
253                 double begClock = clock();
254                 unsigned long int begTime = time(NULL);
255                 
256                 if(pid == 0){
257
258                         for(int i=1;i<ncpus;i++){
259                                 MPI_Send(&abort, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
260                         }
261                         if(abort == 1){ return 0;       }
262
263                         processors = ncpus;
264                         
265                         m->mothurOut("\nGetting preliminary data...\n");
266                         getSingleLookUp();
267                         getJointLookUp();
268                         
269                         vector<string> flowFileVector;
270                         if(flowFilesFileName != "not found"){
271                                 string fName;
272
273                                 ifstream flowFilesFile;
274                                 m->openInputFile(flowFilesFileName, flowFilesFile);
275                                 while(flowFilesFile){
276                                         flowFilesFile >> fName;
277                                         flowFileVector.push_back(fName);
278                                         m->gobble(flowFilesFile);
279                                 }
280                         }
281                         else{
282                                 flowFileVector.push_back(flowFileName);
283                         }
284                         int numFiles = flowFileVector.size();
285
286                         for(int i=1;i<ncpus;i++){
287                                 MPI_Send(&numFiles, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
288                         }
289                         
290                         for(int i=0;i<numFiles;i++){
291                                 flowFileName = flowFileVector[i];
292
293                                 
294                                 
295                                 m->mothurOut("\n>>>>>\tProcessing " + flowFileName + " (file " + toString(i+1) + " of " + toString(numFiles) + ")\t<<<<<\n");
296                                 m->mothurOut("Reading flowgrams...\n");
297                                 
298                                 getFlowData();
299                                 m->mothurOut("Identifying unique flowgrams...\n");
300                                 getUniques();
301
302                                 m->mothurOut("Calculating distances between flowgrams...\n");
303                                 char fileName[1024];
304                                 strcpy(fileName, flowFileName.c_str());
305
306                                 for(int i=1;i<ncpus;i++){
307                                         MPI_Send(&fileName[0], 1024, MPI_CHAR, i, tag, MPI_COMM_WORLD);
308
309                                         MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
310                                         MPI_Send(&numUniques, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
311                                         MPI_Send(&numFlowCells, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
312                                         MPI_Send(&flowDataIntI[0], numSeqs * numFlowCells, MPI_SHORT, i, tag, MPI_COMM_WORLD);
313                                         MPI_Send(&flowDataPrI[0], numSeqs * numFlowCells, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
314                                         MPI_Send(&mapUniqueToSeq[0], numSeqs, MPI_INT, i, tag, MPI_COMM_WORLD);
315                                         MPI_Send(&mapSeqToUnique[0], numSeqs, MPI_INT, i, tag, MPI_COMM_WORLD);
316                                         MPI_Send(&lengths[0], numSeqs, MPI_INT, i, tag, MPI_COMM_WORLD);
317                                         MPI_Send(&jointLookUp[0], NUMBINS * NUMBINS, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
318                                         MPI_Send(&cutoff, 1, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
319                                 }                               
320                                                         
321                                 string distFileName = flowDistMPI(0, int(sqrt(1.0/float(ncpus)) * numUniques));
322                                 
323                                 int done;
324                                 for(int i=1;i<ncpus;i++){
325                                         MPI_Recv(&done, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
326                                         
327                                         m->appendFiles((distFileName + ".temp." + toString(i)), distFileName);
328                                         remove((distFileName + ".temp." + toString(i)).c_str());
329                                 }
330
331                                 string namesFileName = createNamesFile();
332                                 
333                                 m->mothurOut("\nClustering flowgrams...\n");
334                                 string listFileName = cluster(distFileName, namesFileName);
335
336                                 getOTUData(listFileName);
337                                 initPyroCluster();
338
339                                 for(int i=1;i<ncpus;i++){
340                                         MPI_Send(&numOTUs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
341                                         MPI_Send(&singleLookUp[0], singleLookUp.size(), MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
342                                         MPI_Send(&uniqueFlowgrams[0], numFlowCells * numUniques, MPI_SHORT, i, tag, MPI_COMM_WORLD);
343                                         MPI_Send(&sigma, 1, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
344                                 }
345                                 
346                                 
347                                 double maxDelta = 0;
348                                 int iter = 0;
349                                 
350                                 int numOTUsOnCPU = numOTUs / ncpus;
351                                 int numSeqsOnCPU = numSeqs / ncpus;
352                                 m->mothurOut("\nDenoising flowgrams...\n");
353                                 m->mothurOut("iter\tmaxDelta\tnLL\t\tcycletime\n");
354                                 
355                                 while((maxIters == 0 && maxDelta > minDelta) || iter < MIN_ITER || (maxDelta > minDelta && iter < maxIters)){
356
357                                         double cycClock = clock();
358                                         unsigned long int cycTime = time(NULL);
359                                         fill();
360
361                                         int total = singleTau.size();
362                                         for(int i=1;i<ncpus;i++){
363                                                 MPI_Send(&total, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
364                                                 MPI_Send(&change[0], numOTUs, MPI_SHORT, i, tag, MPI_COMM_WORLD);
365                                                 MPI_Send(&centroids[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
366                                                 
367                                                 MPI_Send(&singleTau[0], total, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
368                                                 MPI_Send(&seqNumber[0], total, MPI_INT, i, tag, MPI_COMM_WORLD);
369                                                 MPI_Send(&seqIndex[0], total, MPI_INT, i, tag, MPI_COMM_WORLD);
370                                                 MPI_Send(&nSeqsPerOTU[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
371                                                 MPI_Send(&cumNumSeqs[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
372                                         }
373                                 
374                                         calcCentroidsDriver(0, numOTUsOnCPU);
375                                         
376                                         for(int i=1;i<ncpus;i++){
377                                                 int otuStart = i * numOTUs / ncpus;
378                                                 int otuStop = (i + 1) * numOTUs / ncpus;
379                                                 
380                                                 vector<int> tempCentroids(numOTUs, 0);
381                                                 vector<short> tempChange(numOTUs, 0);
382                                                 
383                                                 MPI_Recv(&tempCentroids[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
384                                                 MPI_Recv(&tempChange[0], numOTUs, MPI_SHORT, i, tag, MPI_COMM_WORLD, &status);
385                                                 
386                                                 for(int j=otuStart;j<otuStop;j++){
387                                                         centroids[j] = tempCentroids[j];
388                                                         change[j] = tempChange[j];
389                                                 }
390                                         }
391                                                                         
392                                         maxDelta = getNewWeights();
393                                         double nLL = getLikelihood();
394                                         checkCentroids();
395                                         
396                                         for(int i=1;i<ncpus;i++){
397                                                 MPI_Send(&centroids[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
398                                                 MPI_Send(&weight[0], numOTUs, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
399                                                 MPI_Send(&change[0], numOTUs, MPI_SHORT, i, tag, MPI_COMM_WORLD);
400                                         }
401                                         
402                                         calcNewDistancesParent(0, numSeqsOnCPU);
403
404                                         total = singleTau.size();
405
406                                         for(int i=1;i<ncpus;i++){
407                                                 int childTotal;
408                                                 int seqStart = i * numSeqs / ncpus;
409                                                 int seqStop = (i + 1) * numSeqs / ncpus;
410                                                 
411                                                 MPI_Recv(&childTotal, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
412
413                                                 vector<int> childSeqIndex(childTotal, 0);
414                                                 vector<double> childSingleTau(childTotal, 0);
415                                                 vector<double> childDist(numSeqs * numOTUs, 0);
416                                                 vector<int> otuIndex(childTotal, 0);
417                                                 
418                                                 MPI_Recv(&childSeqIndex[0], childTotal, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
419                                                 MPI_Recv(&childSingleTau[0], childTotal, MPI_DOUBLE, i, tag, MPI_COMM_WORLD, &status);
420                                                 MPI_Recv(&childDist[0], numOTUs * numSeqs, MPI_DOUBLE, i, tag, MPI_COMM_WORLD, &status);
421                                                 MPI_Recv(&otuIndex[0], childTotal, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
422
423                                                 int oldTotal = total;
424                                                 total += childTotal;
425                                                 singleTau.resize(total, 0);
426                                                 seqIndex.resize(total, 0);
427                                                 seqNumber.resize(total, 0);
428                                                 
429                                                 int childIndex = 0;
430                                                 
431                                                 for(int j=oldTotal;j<total;j++){
432                                                         int otuI = otuIndex[childIndex];
433                                                         int seqI = childSeqIndex[childIndex];
434                                                         
435                                                         singleTau[j] = childSingleTau[childIndex];
436                                                         
437                                                         aaP[otuI][nSeqsPerOTU[otuI]] = j;
438                                                         aaI[otuI][nSeqsPerOTU[otuI]] = seqI;
439                                                         nSeqsPerOTU[otuI]++;
440                                                         childIndex++;
441                                                 }
442                                                 
443                                                 int index = seqStart * numOTUs;
444                                                 for(int j=seqStart;j<seqStop;j++){
445                                                         for(int k=0;k<numOTUs;k++){
446                                                                 dist[index] = childDist[index];
447                                                                 index++;
448                                                         }
449                                                 }                                       
450                                         }
451                                         
452                                         iter++;
453                                         
454                                         m->mothurOut(toString(iter) + '\t' + toString(maxDelta) + '\t' + toString(nLL) + '\t' + toString(time(NULL) - cycTime) + '\t' + toString((clock() - cycClock)/(double)CLOCKS_PER_SEC) + '\n');                  
455
456                                         if((maxIters == 0 && maxDelta > minDelta) || iter < MIN_ITER || (maxDelta > minDelta && iter < maxIters)){
457                                                 int live = 1;
458                                                 for(int i=1;i<ncpus;i++){
459                                                         MPI_Send(&live, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
460                                                 }
461                                         }
462                                         else{
463                                                 int live = 0;
464                                                 for(int i=1;i<ncpus;i++){
465                                                         MPI_Send(&live, 1, MPI_INT, i, tag, MPI_COMM_WORLD); //send kill command
466                                                 }
467                                         }
468                                         
469                                 }       
470                                 
471                                 m->mothurOut("\nFinalizing...\n");
472                                 fill();
473                                 setOTUs();
474                                 vector<int> otuCounts(numOTUs, 0);
475                                 for(int i=0;i<numSeqs;i++)      {       otuCounts[otuData[i]]++;        }
476                                 calcCentroidsDriver(0, numOTUs);
477                                 writeQualities(otuCounts);
478                                 writeSequences(otuCounts);
479                                 writeNames(otuCounts);
480                                 writeClusters(otuCounts);
481                                 writeGroups();
482                                 
483                                 remove(distFileName.c_str());
484                                 remove(namesFileName.c_str());
485                                 remove(listFileName.c_str());
486                                                                  
487                                 m->mothurOut("Total time to process " + toString(flowFileName) + ":\t" + toString(time(NULL) - begTime) + '\t' + toString((clock() - begClock)/(double)CLOCKS_PER_SEC) + '\n');                 
488                         }
489                 }
490                 else{
491                         int abort = 1;
492                         bool live = 1;
493
494                         MPI_Recv(&abort, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
495                         if(abort){      return 0;       }
496
497                         int numFiles;
498                         MPI_Recv(&numFiles, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
499
500                         for(int i=0;i<numFiles;i++){
501                                 //Now into the pyrodist part
502                                 char fileName[1024];
503                                 MPI_Recv(&fileName, 1024, MPI_CHAR, 0, tag, MPI_COMM_WORLD, &status);
504                                 MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
505                                 MPI_Recv(&numUniques, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
506                                 MPI_Recv(&numFlowCells, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
507                                 
508                                 flowDataIntI.resize(numSeqs * numFlowCells);
509                                 flowDataPrI.resize(numSeqs * numFlowCells);
510                                 mapUniqueToSeq.resize(numSeqs);
511                                 mapSeqToUnique.resize(numSeqs);
512                                 lengths.resize(numSeqs);
513                                 jointLookUp.resize(NUMBINS * NUMBINS);
514
515                                 MPI_Recv(&flowDataIntI[0], numSeqs * numFlowCells, MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
516                                 MPI_Recv(&flowDataPrI[0], numSeqs * numFlowCells, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
517                                 MPI_Recv(&mapUniqueToSeq[0], numSeqs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
518                                 MPI_Recv(&mapSeqToUnique[0], numSeqs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
519                                 MPI_Recv(&lengths[0], numSeqs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
520                                 MPI_Recv(&jointLookUp[0], NUMBINS * NUMBINS, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
521                                 MPI_Recv(&cutoff, 1, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
522
523                                 flowFileName = string(fileName);
524                                 int flowDistStart = int(sqrt(float(pid)/float(ncpus)) * numUniques);
525                                 int flowDistEnd = int(sqrt(float(pid+1)/float(ncpus)) * numUniques);
526                                 
527                                 string distanceStringChild = flowDistMPI(flowDistStart, flowDistEnd);
528
529                                 int done = 1;
530                                 MPI_Send(&done, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
531
532                                 //Now into the pyronoise part
533                                 MPI_Recv(&numOTUs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
534                                 
535                                 singleLookUp.resize(HOMOPS * NUMBINS);
536                                 uniqueFlowgrams.resize(numUniques * numFlowCells);
537                                 weight.resize(numOTUs);
538                                 centroids.resize(numOTUs);
539                                 change.resize(numOTUs);
540                                 dist.assign(numOTUs * numSeqs, 0);
541                                 nSeqsPerOTU.resize(numOTUs);
542                                 cumNumSeqs.resize(numOTUs);
543
544                                 MPI_Recv(&singleLookUp[0], singleLookUp.size(), MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
545                                 MPI_Recv(&uniqueFlowgrams[0], uniqueFlowgrams.size(), MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
546                                 MPI_Recv(&sigma, 1, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
547                                 
548                                 int startOTU = pid * numOTUs / ncpus;
549                                 int endOTU = (pid + 1) * numOTUs / ncpus;
550
551                                 int startSeq = pid * numSeqs / ncpus;
552                                 int endSeq = (pid + 1) * numSeqs /ncpus;
553                                 
554                                 int total;
555
556                                 while(live){
557
558                                         MPI_Recv(&total, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
559                                         singleTau.assign(total, 0.0000);
560                                         seqNumber.assign(total, 0);
561                                         seqIndex.assign(total, 0);
562                                         
563                                         MPI_Recv(&change[0], numOTUs, MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
564                                         MPI_Recv(&centroids[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
565                                         MPI_Recv(&singleTau[0], total, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
566                                         MPI_Recv(&seqNumber[0], total, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
567                                         MPI_Recv(&seqIndex[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
568                                         MPI_Recv(&nSeqsPerOTU[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
569                                         MPI_Recv(&cumNumSeqs[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
570                                         
571                                         calcCentroidsDriver(startOTU, endOTU);
572
573                                         MPI_Send(&centroids[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD);
574                                         MPI_Send(&change[0], numOTUs, MPI_SHORT, 0, tag, MPI_COMM_WORLD);
575
576                                         
577                                         MPI_Recv(&centroids[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
578                                         MPI_Recv(&weight[0], numOTUs, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
579                                         MPI_Recv(&change[0], numOTUs, MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
580
581                                         vector<int> otuIndex(total, 0);
582                                         calcNewDistancesChildMPI(startSeq, endSeq, otuIndex);
583                                         total = otuIndex.size();
584                                         
585                                         MPI_Send(&total, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
586                                         MPI_Send(&seqIndex[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD);
587                                         MPI_Send(&singleTau[0], total, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD);
588                                         MPI_Send(&dist[0], numOTUs * numSeqs, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD);
589                                         MPI_Send(&otuIndex[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD);
590                                 
591                                         MPI_Recv(&live, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
592                                 }
593                         }
594                 }               
595
596                 MPI_Barrier(MPI_COMM_WORLD);
597                 return 0;
598
599         }
600         catch(exception& e) {
601                 m->errorOut(e, "ShhherCommand", "execute");
602                 exit(1);
603         }
604 }
605
606 /**************************************************************************************************/
607
608 string ShhherCommand::flowDistMPI(int startSeq, int stopSeq){
609         try{            
610                 ostringstream outStream;
611                 outStream.setf(ios::fixed, ios::floatfield);
612                 outStream.setf(ios::dec, ios::basefield);
613                 outStream.setf(ios::showpoint);
614                 outStream.precision(6);
615                 
616                 int begTime = time(NULL);
617                 double begClock = clock();
618                 
619                 for(int i=startSeq;i<stopSeq;i++){
620                         for(int j=0;j<i;j++){
621                                 float flowDistance = calcPairwiseDist(mapUniqueToSeq[i], mapUniqueToSeq[j]);
622                                 
623                                 if(flowDistance < 1e-6){
624                                         outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << 0.000000 << endl;
625                                 }
626                                 else if(flowDistance <= cutoff){
627                                         outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << flowDistance << endl;
628                                 }
629                         }
630                         if(i % 100 == 0){
631                                 m->mothurOut(toString(i) + '\t' + toString(time(NULL) - begTime) + '\t' + toString((clock()-begClock)/CLOCKS_PER_SEC) + '\n');
632                         }
633                 }
634                 
635                 m->mothurOut(toString(stopSeq) + '\t' + toString(time(NULL) - begTime) + '\t' + toString((clock()-begClock)/CLOCKS_PER_SEC) + '\n');
636                 
637                 string fDistFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.dist";
638                 if(pid != 0){   fDistFileName += ".temp." + toString(pid);      }
639
640                 ofstream distFile(fDistFileName.c_str());
641                 distFile << outStream.str();            
642                 distFile.close();
643                 
644                 return fDistFileName;
645         }
646         catch(exception& e) {
647                 m->errorOut(e, "ShhherCommand", "flowDistParentFork");
648                 exit(1);
649         }
650 }
651
652 #else
653 //**********************************************************************************************************************
654
655 int ShhherCommand::execute(){
656         try {
657                 if (abort == true) { return 0; }
658                 
659                 getSingleLookUp();
660                 getJointLookUp();
661                                 
662                 vector<string> flowFileVector;
663                 if(flowFilesFileName != "not found"){
664                         string fName;
665                         
666                         ifstream flowFilesFile;
667                         m->openInputFile(flowFilesFileName, flowFilesFile);
668                         while(flowFilesFile){
669                                 flowFilesFile >> fName;
670                                 flowFileVector.push_back(fName);
671                                 m->gobble(flowFilesFile);
672                         }
673                 }
674                 else{
675                         flowFileVector.push_back(flowFileName);
676                 }
677                 int numFiles = flowFileVector.size();
678                 
679                 
680                 for(int i=0;i<numFiles;i++){
681                         flowFileName = flowFileVector[i];
682
683                         m->mothurOut("\n>>>>>\tProcessing " + flowFileName + " (file " + toString(i+1) + " of " + toString(numFiles) + ")\t<<<<<\n");
684                         m->mothurOut("Reading flowgrams...\n");
685                         getFlowData();
686                         
687                         m->mothurOut("Identifying unique flowgrams...\n");
688                         getUniques();
689                         
690                         
691                         m->mothurOut("Calculating distances between flowgrams...\n");
692                         string distFileName = createDistFile(processors);
693                         string namesFileName = createNamesFile();
694                                 
695                         m->mothurOut("\nClustering flowgrams...\n");
696                         string listFileName = cluster(distFileName, namesFileName);
697                         getOTUData(listFileName);
698                         
699                         initPyroCluster();
700                         
701                         double maxDelta = 0;
702                         int iter = 0;
703                         
704                         double begClock = clock();
705                         unsigned long int begTime = time(NULL);
706
707                         m->mothurOut("\nDenoising flowgrams...\n");
708                         m->mothurOut("iter\tmaxDelta\tnLL\t\tcycletime\n");
709                         
710                         while((maxIters == 0 && maxDelta > minDelta) || iter < MIN_ITER || (maxDelta > minDelta && iter < maxIters)){
711                                 
712                                 double cycClock = clock();
713                                 unsigned long int cycTime = time(NULL);
714                                 fill();
715                                 
716                                 calcCentroids();
717                                 
718                                 maxDelta = getNewWeights();
719                                 double nLL = getLikelihood();
720                                 checkCentroids();
721                                 
722                                 calcNewDistances();
723
724                                 iter++;
725                                 
726                                 m->mothurOut(toString(iter) + '\t' + toString(maxDelta) + '\t' + toString(nLL) + '\t' + toString(time(NULL) - cycTime) + '\t' + toString((clock() - cycClock)/(double)CLOCKS_PER_SEC) + '\n');
727
728                         }       
729                         
730                         m->mothurOut("\nFinalizing...\n");
731                         fill();
732                         setOTUs();
733                         
734                         vector<int> otuCounts(numOTUs, 0);
735                         for(int i=0;i<numSeqs;i++)      {       otuCounts[otuData[i]]++;        }
736                         
737                         calcCentroidsDriver(0, numOTUs);
738                         writeQualities(otuCounts);
739                         writeSequences(otuCounts);
740                         writeNames(otuCounts);
741                         writeClusters(otuCounts);
742                         writeGroups();
743                         
744                         remove(distFileName.c_str());
745                         remove(namesFileName.c_str());
746                         remove(listFileName.c_str());
747                         
748                         m->mothurOut("Total time to process " + flowFileName + ":\t" + toString(time(NULL) - begTime) + '\t' + toString((clock() - begClock)/(double)CLOCKS_PER_SEC) + '\n');
749                 }
750                 return 0;
751         }
752         catch(exception& e) {
753                 m->errorOut(e, "ShhherCommand", "execute");
754                 exit(1);
755         }
756 }
757 #endif
758 /**************************************************************************************************/
759
760 void ShhherCommand::getFlowData(){
761         try{
762                 ifstream flowFile;
763                 m->openInputFile(flowFileName, flowFile);
764                 
765                 string seqName;
766                 seqNameVector.clear();
767                 lengths.clear();
768                 flowDataIntI.clear();
769                 nameMap.clear();
770                 
771                 
772                 int currentNumFlowCells;
773                 
774                 float intensity;
775                 
776                 flowFile >> numFlowCells;
777                 int index = 0;//pcluster
778                 while(!flowFile.eof()){
779                         flowFile >> seqName >> currentNumFlowCells;
780                         lengths.push_back(currentNumFlowCells);
781
782                         seqNameVector.push_back(seqName);
783                         nameMap[seqName] = index++;//pcluster
784
785                         for(int i=0;i<numFlowCells;i++){
786                                 flowFile >> intensity;
787                                 if(intensity > 9.99)    {       intensity = 9.99;       }
788                                 int intI = int(100 * intensity + 0.0001);
789                                 flowDataIntI.push_back(intI);
790                         }
791                         m->gobble(flowFile);
792                 }
793                 flowFile.close();
794                 
795                 numSeqs = seqNameVector.size();         
796                 
797                 for(int i=0;i<numSeqs;i++){
798                         int iNumFlowCells = i * numFlowCells;
799                         for(int j=lengths[i];j<numFlowCells;j++){
800                                 flowDataIntI[iNumFlowCells + j] = 0;
801                         }
802                 }
803                 
804         }
805         catch(exception& e) {
806                 m->errorOut(e, "ShhherCommand", "getFlowData");
807                 exit(1);
808         }
809 }
810
811 /**************************************************************************************************/
812
813 void ShhherCommand::getSingleLookUp(){
814         try{
815                 //      these are the -log probabilities that a signal corresponds to a particular homopolymer length
816                 singleLookUp.assign(HOMOPS * NUMBINS, 0);
817                 
818                 int index = 0;
819                 ifstream lookUpFile;
820                 m->openInputFile(lookupFileName, lookUpFile);
821                 
822                 for(int i=0;i<HOMOPS;i++){
823                         float logFracFreq;
824                         lookUpFile >> logFracFreq;
825                         
826                         for(int j=0;j<NUMBINS;j++)      {
827                                 lookUpFile >> singleLookUp[index];
828                                 index++;
829                         }
830                 }       
831                 lookUpFile.close();
832         }
833         catch(exception& e) {
834                 m->errorOut(e, "ShhherCommand", "getSingleLookUp");
835                 exit(1);
836         }
837 }
838
839 /**************************************************************************************************/
840
841 void ShhherCommand::getJointLookUp(){
842         try{
843                 
844                 //      the most likely joint probability (-log) that two intenities have the same polymer length
845                 jointLookUp.resize(NUMBINS * NUMBINS, 0);
846                 
847                 for(int i=0;i<NUMBINS;i++){
848                         for(int j=0;j<NUMBINS;j++){             
849                                 
850                                 double minSum = 100000000;
851                                 
852                                 for(int k=0;k<HOMOPS;k++){
853                                         double sum = singleLookUp[k * NUMBINS + i] + singleLookUp[k * NUMBINS + j];
854                                         
855                                         if(sum < minSum)        {       minSum = sum;           }
856                                 }       
857                                 jointLookUp[i * NUMBINS + j] = minSum;
858                         }
859                 }
860         }
861         catch(exception& e) {
862                 m->errorOut(e, "ShhherCommand", "getJointLookUp");
863                 exit(1);
864         }
865 }
866
867 /**************************************************************************************************/
868
869 double ShhherCommand::getProbIntensity(int intIntensity){                          
870         try{
871
872                 double minNegLogProb = 100000000; 
873
874                 
875                 for(int i=0;i<HOMOPS;i++){//loop signal strength
876                         float negLogProb = singleLookUp[i * NUMBINS + intIntensity];
877                         if(negLogProb < minNegLogProb)  {       minNegLogProb = negLogProb; }
878                 }
879                 
880                 return minNegLogProb;
881         }
882         catch(exception& e) {
883                 m->errorOut(e, "ShhherCommand", "getProbIntensity");
884                 exit(1);
885         }
886 }
887
888 /**************************************************************************************************/
889
890 void ShhherCommand::getUniques(){
891         try{
892                 
893                 
894                 numUniques = 0;
895                 uniqueFlowgrams.assign(numFlowCells * numSeqs, -1);
896                 uniqueCount.assign(numSeqs, 0);                                                 //      anWeights
897                 uniqueLengths.assign(numSeqs, 0);
898                 mapSeqToUnique.assign(numSeqs, -1);
899                 mapUniqueToSeq.assign(numSeqs, -1);
900                 
901                 vector<short> uniqueFlowDataIntI(numFlowCells * numSeqs, -1);
902                 
903                 for(int i=0;i<numSeqs;i++){
904                         int index = 0;
905                         
906                         vector<short> current(numFlowCells);
907                         for(int j=0;j<numFlowCells;j++){
908                                 current[j] = short(((flowDataIntI[i * numFlowCells + j] + 50.0)/100.0));
909                         }
910                                                 
911                         for(int j=0;j<numUniques;j++){
912                                 int offset = j * numFlowCells;
913                                 bool toEnd = 1;
914                                 
915                                 int shorterLength;
916                                 if(lengths[i] < uniqueLengths[j])       {       shorterLength = lengths[i];                     }
917                                 else                                                            {       shorterLength = uniqueLengths[j];       }
918
919                                 for(int k=0;k<shorterLength;k++){
920                                         if(current[k] != uniqueFlowgrams[offset + k]){
921                                                 toEnd = 0;
922                                                 break;
923                                         }
924                                 }
925                                 
926                                 if(toEnd){
927                                         mapSeqToUnique[i] = j;
928                                         uniqueCount[j]++;
929                                         index = j;
930                                         if(lengths[i] > uniqueLengths[j])       {       uniqueLengths[j] = lengths[i];  }
931                                         break;
932                                 }
933                                 index++;
934                         }
935                         
936                         if(index == numUniques){
937                                 uniqueLengths[numUniques] = lengths[i];
938                                 uniqueCount[numUniques] = 1;
939                                 mapSeqToUnique[i] = numUniques;//anMap
940                                 mapUniqueToSeq[numUniques] = i;//anF
941                                 
942                                 for(int k=0;k<numFlowCells;k++){
943                                         uniqueFlowgrams[numUniques * numFlowCells + k] = current[k];
944                                         uniqueFlowDataIntI[numUniques * numFlowCells + k] = flowDataIntI[i * numFlowCells + k];
945                                 }
946                                 
947                                 numUniques++;
948                         }
949                 }
950                 uniqueFlowDataIntI.resize(numFlowCells * numUniques);
951                 uniqueLengths.resize(numUniques);       
952                 
953                 flowDataPrI.resize(numSeqs * numFlowCells, 0);
954                 for(int i=0;i<flowDataPrI.size();i++)   {       flowDataPrI[i] = getProbIntensity(flowDataIntI[i]);             }
955         }
956         catch(exception& e) {
957                 m->errorOut(e, "ShhherCommand", "getUniques");
958                 exit(1);
959         }
960 }
961
962 /**************************************************************************************************/
963
964 float ShhherCommand::calcPairwiseDist(int seqA, int seqB){
965         try{
966                 int minLength = lengths[mapSeqToUnique[seqA]];
967                 if(lengths[seqB] < minLength){  minLength = lengths[mapSeqToUnique[seqB]];      }
968                 
969                 int ANumFlowCells = seqA * numFlowCells;
970                 int BNumFlowCells = seqB * numFlowCells;
971                 
972                 float dist = 0;
973                 
974                 for(int i=0;i<minLength;i++){
975                         int flowAIntI = flowDataIntI[ANumFlowCells + i];
976                         float flowAPrI = flowDataPrI[ANumFlowCells + i];
977                         
978                         int flowBIntI = flowDataIntI[BNumFlowCells + i];
979                         float flowBPrI = flowDataPrI[BNumFlowCells + i];
980                         dist += jointLookUp[flowAIntI * NUMBINS + flowBIntI] - flowAPrI - flowBPrI;
981                 }
982                 
983                 dist /= (float) minLength;
984                 return dist;
985         }
986         catch(exception& e) {
987                 m->errorOut(e, "ShhherCommand", "calcPairwiseDist");
988                 exit(1);
989         }
990 }
991
992 /**************************************************************************************************/
993
994 void ShhherCommand::flowDistParentFork(string distFileName, int startSeq, int stopSeq){
995         try{            
996
997                 ostringstream outStream;
998                 outStream.setf(ios::fixed, ios::floatfield);
999                 outStream.setf(ios::dec, ios::basefield);
1000                 outStream.setf(ios::showpoint);
1001                 outStream.precision(6);
1002                 
1003                 int begTime = time(NULL);
1004                 double begClock = clock();
1005
1006                 for(int i=startSeq;i<stopSeq;i++){
1007                         for(int j=0;j<i;j++){
1008                                 float flowDistance = calcPairwiseDist(mapUniqueToSeq[i], mapUniqueToSeq[j]);
1009
1010                                 if(flowDistance < 1e-6){
1011                                         outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << 0.000000 << endl;
1012                                 }
1013                                 else if(flowDistance <= cutoff){
1014                                         outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << flowDistance << endl;
1015                                 }
1016                         }
1017                         if(i % 100 == 0){
1018                                 m->mothurOut(toString(i) + "\t" + toString(time(NULL) - begTime));
1019                                 m->mothurOut("\t" + toString((clock()-begClock)/CLOCKS_PER_SEC));
1020                                 m->mothurOutEndLine();
1021                         }
1022                 }
1023                 m->mothurOut(toString(stopSeq-1) + "\t" + toString(time(NULL) - begTime));
1024                 m->mothurOut("\t" + toString((clock()-begClock)/CLOCKS_PER_SEC));
1025                 m->mothurOutEndLine();
1026                 
1027                 ofstream distFile(distFileName.c_str());
1028                 distFile << outStream.str();            
1029                 distFile.close();
1030         }
1031         catch(exception& e) {
1032                 m->errorOut(e, "ShhherCommand", "flowDistParentFork");
1033                 exit(1);
1034         }
1035 }
1036
1037 /**************************************************************************************************/
1038
1039 string ShhherCommand::createDistFile(int processors){
1040         try{
1041                 string fDistFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.dist";
1042                                 
1043                 unsigned long int begTime = time(NULL);
1044                 double begClock = clock();
1045
1046                 vector<int> start;
1047                 vector<int> end;
1048                 
1049 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1050                 if(processors == 1)     {       flowDistParentFork(fDistFileName, 0, numUniques);               }
1051                 else{ //you have multiple processors
1052                         
1053                         if (numSeqs < processors){      processors = 1; }
1054                         
1055                         vector<int> start(processors, 0);
1056                         vector<int> end(processors, 0);
1057                         
1058                         for (int i = 0; i < processors; i++) {
1059                                 start[i] = int(sqrt(float(i)/float(processors)) * numUniques);
1060                                 end[i] = int(sqrt(float(i+1)/float(processors)) * numUniques);
1061                         }
1062                         
1063                         int process = 1;
1064                         vector<int> processIDs;
1065                         
1066                         //loop through and create all the processes you want
1067                         while (process != processors) {
1068                                 int pid = fork();
1069                                 
1070                                 if (pid > 0) {
1071                                         processIDs.push_back(pid);  //create map from line number to pid so you can append files in correct order later
1072                                         process++;
1073                                 }else if (pid == 0){
1074                                         flowDistParentFork(fDistFileName + toString(getpid()) + ".temp", start[process], end[process]);
1075                                         exit(0);
1076                                 }else { 
1077                                         m->mothurOut("[ERROR]: unable to spawn the necessary processes. Error code: " + toString(pid)); m->mothurOutEndLine(); 
1078                                         perror(" : ");
1079                                         for (int i=0;i<processIDs.size();i++) {  int temp = processIDs[i]; kill (temp, SIGINT); }
1080                                         exit(0);
1081                                 }
1082                         }
1083                         
1084                         //parent does its part
1085                         flowDistParentFork(fDistFileName, start[0], end[0]);
1086                         
1087                         //force parent to wait until all the processes are done
1088                         for (int i=0;i<processIDs.size();i++) { 
1089                                 int temp = processIDs[i];
1090                                 wait(&temp);
1091                         }
1092                         
1093                         //append and remove temp files
1094                         for (int i=0;i<processIDs.size();i++) { 
1095                                 m->appendFiles((fDistFileName + toString(processIDs[i]) + ".temp"), fDistFileName);
1096                                 remove((fDistFileName + toString(processIDs[i]) + ".temp").c_str());
1097                         }
1098                         
1099                 }
1100                 
1101 #else
1102                 flowDistParentFork(fDistFileName, 0, numUniques);
1103 #endif
1104
1105                 m->mothurOutEndLine();
1106                 
1107                 m->mothurOut("Total time: " + toString(time(NULL) - begTime) + '\t' + toString((clock() - begClock)/CLOCKS_PER_SEC) + '\n');
1108                 
1109
1110                 return fDistFileName;
1111         }
1112         catch(exception& e) {
1113                 m->errorOut(e, "ShhherCommand", "createDistFile");
1114                 exit(1);
1115         }
1116         
1117 }
1118
1119 /**************************************************************************************************/
1120
1121 string ShhherCommand::createNamesFile(){
1122         try{
1123                 
1124                 vector<string> duplicateNames(numUniques, "");
1125                 for(int i=0;i<numSeqs;i++){
1126                         duplicateNames[mapSeqToUnique[i]] += seqNameVector[i] + ',';
1127                 }
1128                 
1129                 string nameFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.names";
1130                 
1131                 ofstream nameFile;
1132                 m->openOutputFile(nameFileName, nameFile);
1133                 
1134                 for(int i=0;i<numUniques;i++){
1135 //                      nameFile << seqNameVector[mapUniqueToSeq[i]] << '\t' << duplicateNames[i].substr(0, duplicateNames[i].find_last_of(',')) << endl;
1136                         nameFile << mapUniqueToSeq[i] << '\t' << duplicateNames[i].substr(0, duplicateNames[i].find_last_of(',')) << endl;
1137                 }
1138                 
1139                 nameFile.close();
1140                 return  nameFileName;
1141         }
1142         catch(exception& e) {
1143                 m->errorOut(e, "ShhherCommand", "createNamesFile");
1144                 exit(1);
1145         }
1146 }
1147
1148 //**********************************************************************************************************************
1149
1150 string ShhherCommand::cluster(string distFileName, string namesFileName){
1151         try {
1152                 
1153                 
1154                 globaldata->setNameFile(namesFileName);
1155                 globaldata->setColumnFile(distFileName);
1156                 globaldata->setFormat("column");
1157                 
1158                 ReadMatrix* read = new ReadColumnMatrix(distFileName);  
1159                 read->setCutoff(cutoff);
1160                 
1161                 NameAssignment* clusterNameMap = new NameAssignment(namesFileName);
1162                 clusterNameMap->readMap();
1163                 read->read(clusterNameMap);
1164                 
1165                 ListVector* list = read->getListVector();
1166                 SparseMatrix* matrix = read->getMatrix();
1167                 
1168                 delete read; 
1169                 delete clusterNameMap; 
1170                                 
1171                 RAbundVector* rabund = new RAbundVector(list->getRAbundVector());
1172                 
1173                 Cluster* cluster = new CompleteLinkage(rabund, list, matrix, cutoff, "furthest"); 
1174                 string tag = cluster->getTag();
1175                 
1176                 double clusterCutoff = cutoff;
1177                 while (matrix->getSmallDist() <= clusterCutoff && matrix->getNNodes() > 0){
1178                         cluster->update(clusterCutoff);
1179                 }
1180                 
1181                 list->setLabel(toString(cutoff));
1182                 
1183                 string listFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.list";
1184                 ofstream listFile;
1185                 m->openOutputFile(listFileName, listFile);
1186                 list->print(listFile);
1187                 listFile.close();
1188                 
1189                 delete matrix;  delete cluster; delete rabund; delete list;
1190         
1191                 return listFileName;
1192         }
1193         catch(exception& e) {
1194                 m->errorOut(e, "ShhherCommand", "cluster");
1195                 exit(1);        
1196         }               
1197 }
1198
1199 /**************************************************************************************************/
1200
1201 void ShhherCommand::getOTUData(string listFileName){
1202         try {
1203
1204                 ifstream listFile;
1205                 m->openInputFile(listFileName, listFile);
1206                 string label;
1207                 
1208                 listFile >> label >> numOTUs;
1209
1210                 otuData.assign(numSeqs, 0);
1211                 cumNumSeqs.assign(numOTUs, 0);
1212                 nSeqsPerOTU.assign(numOTUs, 0);
1213                 aaP.resize(numOTUs);
1214                 
1215                 seqNumber.clear();
1216                 aaI.clear();
1217                 seqIndex.clear();
1218                 
1219                 string singleOTU = "";
1220                 
1221                 for(int i=0;i<numOTUs;i++){
1222
1223                         listFile >> singleOTU;
1224                         
1225                         istringstream otuString(singleOTU);
1226
1227                         while(otuString){
1228                                 
1229                                 string seqName = "";
1230                                 
1231                                 for(int j=0;j<singleOTU.length();j++){
1232                                         char letter = otuString.get();
1233                                         
1234                                         if(letter != ','){
1235                                                 seqName += letter;
1236                                         }
1237                                         else{
1238                                                 map<string,int>::iterator nmIt = nameMap.find(seqName);
1239                                                 int index = nmIt->second;
1240                                                 
1241                                                 nameMap.erase(nmIt);
1242                                                 
1243                                                 otuData[index] = i;
1244                                                 nSeqsPerOTU[i]++;
1245                                                 aaP[i].push_back(index);
1246                                                 seqName = "";
1247                                         }
1248                                 }
1249                                 
1250                                 map<string,int>::iterator nmIt = nameMap.find(seqName);
1251
1252                                 int index = nmIt->second;
1253                                 nameMap.erase(nmIt);
1254
1255                                 otuData[index] = i;
1256                                 nSeqsPerOTU[i]++;
1257                                 aaP[i].push_back(index);        
1258                                 
1259                                 otuString.get();
1260                         }
1261                         
1262                         sort(aaP[i].begin(), aaP[i].end());
1263                         for(int j=0;j<nSeqsPerOTU[i];j++){
1264                                 seqNumber.push_back(aaP[i][j]);
1265                         }
1266                         for(int j=nSeqsPerOTU[i];j<numSeqs;j++){
1267                                 aaP[i].push_back(0);
1268                         }
1269                 }
1270                 
1271                 for(int i=1;i<numOTUs;i++){
1272                         cumNumSeqs[i] = cumNumSeqs[i-1] + nSeqsPerOTU[i-1];
1273                 }
1274                 aaI = aaP;
1275                 seqIndex = seqNumber;
1276                 
1277                 listFile.close();       
1278         }
1279         catch(exception& e) {
1280                 m->errorOut(e, "ShhherCommand", "getOTUData");
1281                 exit(1);        
1282         }               
1283 }
1284
1285 /**************************************************************************************************/
1286
1287 void ShhherCommand::initPyroCluster(){                          
1288         try{
1289                 dist.assign(numSeqs * numOTUs, 0);
1290                 change.assign(numOTUs, 1);
1291                 centroids.assign(numOTUs, -1);
1292                 weight.assign(numOTUs, 0);
1293                 singleTau.assign(numSeqs, 1.0);
1294                 
1295                 nSeqsBreaks.assign(processors+1, 0);
1296                 nOTUsBreaks.assign(processors+1, 0);
1297
1298                 nSeqsBreaks[0] = 0;
1299                 for(int i=0;i<processors;i++){
1300                         nSeqsBreaks[i+1] = nSeqsBreaks[i] + (int)((double) numSeqs / (double) processors);
1301                         nOTUsBreaks[i+1] = nOTUsBreaks[i] + (int)((double) numOTUs / (double) processors);
1302                 }
1303                 nSeqsBreaks[processors] = numSeqs;
1304                 nOTUsBreaks[processors] = numOTUs;
1305         }
1306         catch(exception& e) {
1307                 m->errorOut(e, "ShhherCommand", "initPyroCluster");
1308                 exit(1);        
1309         }
1310 }
1311
1312 /**************************************************************************************************/
1313
1314 void ShhherCommand::fill(){
1315         try {
1316                 int index = 0;
1317                 for(int i=0;i<numOTUs;i++){
1318                         cumNumSeqs[i] = index;
1319                         for(int j=0;j<nSeqsPerOTU[i];j++){
1320                                 seqNumber[index] = aaP[i][j];
1321                                 seqIndex[index] = aaI[i][j];
1322                                 
1323                                 index++;
1324                         }
1325                 }
1326         }
1327         catch(exception& e) {
1328                 m->errorOut(e, "ShhherCommand", "fill");
1329                 exit(1);        
1330         }               
1331 }
1332
1333 /**************************************************************************************************/
1334
1335 void ShhherCommand::calcCentroids(){                          
1336         try{
1337                 
1338 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1339
1340                 if(processors == 1)     {
1341                         calcCentroidsDriver(0, numOTUs);                
1342                 }
1343                 else{ //you have multiple processors
1344                         if (numOTUs < processors){      processors = 1; }
1345                         
1346                         int process = 1;
1347                         vector<int> processIDs;
1348                         
1349                         //loop through and create all the processes you want
1350                         while (process != processors) {
1351                                 int pid = vfork();
1352                                 
1353                                 if (pid > 0) {
1354                                         processIDs.push_back(pid);  //create map from line number to pid so you can append files in correct order later
1355                                         process++;
1356                                 }else if (pid == 0){
1357                                         calcCentroidsDriver(nOTUsBreaks[process], nOTUsBreaks[process+1]);
1358                                         exit(0);
1359                                 }else { 
1360                                         m->mothurOut("[ERROR]: unable to spawn the necessary processes. Error code: " + toString(pid)); m->mothurOutEndLine(); 
1361                                         perror(" : ");
1362                                         for (int i=0;i<processIDs.size();i++) {  int temp = processIDs[i]; kill (temp, SIGINT); }
1363                                         exit(0);
1364                                 }
1365                         }
1366                         
1367                         //parent does its part
1368                         calcCentroidsDriver(nOTUsBreaks[0], nOTUsBreaks[1]);
1369
1370                         //force parent to wait until all the processes are done
1371                         for (int i=0;i<processIDs.size();i++) { 
1372                                 int temp = processIDs[i];
1373                                 wait(&temp);
1374                         }
1375                 }
1376                 
1377 #else
1378                 calcCentroidsDriver(0, numOTUs);
1379 #endif          
1380         }
1381         catch(exception& e) {
1382                 m->errorOut(e, "ShhherCommand", "calcCentroidsDriver");
1383                 exit(1);        
1384         }               
1385 }
1386
1387 /**************************************************************************************************/
1388
1389 void ShhherCommand::calcCentroidsDriver(int start, int finish){                          
1390         
1391         //this function gets the most likely homopolymer length at a flow position for a group of sequences
1392         //within an otu
1393         
1394         try{
1395                 
1396                 for(int i=start;i<finish;i++){
1397                         
1398                         double count = 0;
1399                         int position = 0;
1400                         int minFlowGram = 100000000;
1401                         double minFlowValue = 1e8;
1402                         change[i] = 0; //FALSE
1403                         
1404                         for(int j=0;j<nSeqsPerOTU[i];j++){
1405                                 count += singleTau[seqNumber[cumNumSeqs[i] + j]];
1406                         }
1407                         
1408                         if(nSeqsPerOTU[i] > 0 && count > MIN_COUNT){
1409                                 vector<double> adF(nSeqsPerOTU[i]);
1410                                 vector<int> anL(nSeqsPerOTU[i]);
1411                                 
1412                                 for(int j=0;j<nSeqsPerOTU[i];j++){
1413                                         int index = cumNumSeqs[i] + j;
1414                                         int nI = seqIndex[index];
1415                                         int nIU = mapSeqToUnique[nI];
1416                                         
1417                                         int k;
1418                                         for(k=0;k<position;k++){
1419                                                 if(nIU == anL[k]){
1420                                                         break;
1421                                                 }
1422                                         }
1423                                         if(k == position){
1424                                                 anL[position] = nIU;
1425                                                 adF[position] = 0.0000;
1426                                                 position++;
1427                                         }                                               
1428                                 }
1429                                 
1430                                 for(int j=0;j<nSeqsPerOTU[i];j++){
1431                                         int index = cumNumSeqs[i] + j;
1432                                         int nI = seqIndex[index];
1433                                         
1434                                         double tauValue = singleTau[seqNumber[index]];
1435                                         
1436                                         for(int k=0;k<position;k++){
1437                                                 double dist = getDistToCentroid(anL[k], nI, lengths[nI]);
1438                                                 adF[k] += dist * tauValue;
1439                                         }
1440                                 }
1441                                 
1442                                 for(int j=0;j<position;j++){
1443                                         if(adF[j] < minFlowValue){
1444                                                 minFlowGram = j;
1445                                                 minFlowValue = adF[j];
1446                                         }
1447                                 }
1448                                 
1449                                 if(centroids[i] != anL[minFlowGram]){
1450                                         change[i] = 1;
1451                                         centroids[i] = anL[minFlowGram];
1452                                 }
1453                         }
1454                         else if(centroids[i] != -1){
1455                                 change[i] = 1;
1456                                 centroids[i] = -1;                      
1457                         }
1458                 }
1459         }
1460         catch(exception& e) {
1461                 m->errorOut(e, "ShhherCommand", "calcCentroidsDriver");
1462                 exit(1);        
1463         }               
1464 }
1465
1466 /**************************************************************************************************/
1467
1468 double ShhherCommand::getDistToCentroid(int cent, int flow, int length){
1469         try{
1470                 
1471                 int flowAValue = cent * numFlowCells;
1472                 int flowBValue = flow * numFlowCells;
1473                 
1474                 double dist = 0;
1475                 
1476                 for(int i=0;i<length;i++){
1477                         dist += singleLookUp[uniqueFlowgrams[flowAValue] * NUMBINS + flowDataIntI[flowBValue]];
1478                         flowAValue++;
1479                         flowBValue++;
1480                 }
1481                 
1482                 return dist / (double)length;
1483         }
1484         catch(exception& e) {
1485                 m->errorOut(e, "ShhherCommand", "getDistToCentroid");
1486                 exit(1);        
1487         }               
1488 }
1489
1490 /**************************************************************************************************/
1491
1492 double ShhherCommand::getNewWeights(){
1493         try{
1494                 
1495                 double maxChange = 0;
1496                 
1497                 for(int i=0;i<numOTUs;i++){
1498                         
1499                         double difference = weight[i];
1500                         weight[i] = 0;
1501                         
1502                         for(int j=0;j<nSeqsPerOTU[i];j++){
1503                                 int index = cumNumSeqs[i] + j;
1504                                 double tauValue = singleTau[seqNumber[index]];
1505                                 weight[i] += tauValue;
1506                         }
1507                         
1508                         difference = fabs(weight[i] - difference);
1509                         if(difference > maxChange){     maxChange = difference; }
1510                 }
1511                 return maxChange;
1512         }
1513         catch(exception& e) {
1514                 m->errorOut(e, "ShhherCommand", "getNewWeights");
1515                 exit(1);        
1516         }               
1517 }
1518
1519 /**************************************************************************************************/
1520
1521 double ShhherCommand::getLikelihood(){
1522         
1523         try{
1524                 
1525                 vector<long double> P(numSeqs, 0);
1526                 int effNumOTUs = 0;
1527                 
1528                 for(int i=0;i<numOTUs;i++){
1529                         if(weight[i] > MIN_WEIGHT){
1530                                 effNumOTUs++;
1531                         }
1532                 }
1533                 
1534                 string hold;
1535                 for(int i=0;i<numOTUs;i++){
1536                         for(int j=0;j<nSeqsPerOTU[i];j++){
1537                                 int index = cumNumSeqs[i] + j;
1538                                 int nI = seqIndex[index];
1539                                 double singleDist = dist[seqNumber[index]];
1540                                 
1541                                 P[nI] += weight[i] * exp(-singleDist * sigma);
1542                         }
1543                 }
1544                 double nLL = 0.00;
1545                 for(int i=0;i<numSeqs;i++){
1546                         if(P[i] == 0){  P[i] = DBL_EPSILON;     }
1547
1548                         nLL += -log(P[i]);
1549                 }
1550                 
1551                 nLL = nLL -(double)numSeqs * log(sigma);
1552
1553                 return nLL; 
1554         }
1555         catch(exception& e) {
1556                 m->errorOut(e, "ShhherCommand", "getNewWeights");
1557                 exit(1);        
1558         }               
1559 }
1560
1561 /**************************************************************************************************/
1562
1563 void ShhherCommand::checkCentroids(){
1564         try{
1565                 vector<int> unique(numOTUs, 1);
1566                 
1567                 for(int i=0;i<numOTUs;i++){
1568                         if(centroids[i] == -1 || weight[i] < MIN_WEIGHT){
1569                                 unique[i] = -1;
1570                         }
1571                 }
1572                 
1573                 for(int i=0;i<numOTUs;i++){
1574                         if(unique[i] == 1){
1575                                 for(int j=i+1;j<numOTUs;j++){
1576                                         if(unique[j] == 1){
1577                                                 
1578                                                 if(centroids[j] == centroids[i]){
1579                                                         unique[j] = 0;
1580                                                         centroids[j] = -1;
1581                                                         
1582                                                         weight[i] += weight[j];
1583                                                         weight[j] = 0.0;
1584                                                 }
1585                                         }
1586                                 }
1587                         }
1588                 }
1589         }
1590         catch(exception& e) {
1591                 m->errorOut(e, "ShhherCommand", "checkCentroids");
1592                 exit(1);        
1593         }               
1594 }
1595
1596 /**************************************************************************************************/
1597
1598 void ShhherCommand::calcNewDistances(){                          
1599         try{
1600                 
1601 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1602                 
1603                 if(processors == 1)     {
1604                         calcNewDistancesParent(0, numSeqs);             
1605                 }
1606                 else{ //you have multiple processors
1607                         if (numSeqs < processors){      processors = 1; }
1608                         
1609                         vector<vector<int> > child_otuIndex(processors);
1610                         vector<vector<int> > child_seqIndex(processors);
1611                         vector<vector<double> > child_singleTau(processors);                    
1612                         vector<int> totals(processors);
1613                         
1614                         int process = 1;
1615                         vector<int> processIDs;
1616
1617                         //loop through and create all the processes you want
1618                         while (process != processors) {
1619                                 int pid = vfork();
1620                                 
1621                                 if (pid > 0) {
1622                                         processIDs.push_back(pid);  //create map from line number to pid so you can append files in correct order later
1623                                         process++;
1624                                 }else if (pid == 0){
1625                                         calcNewDistancesChild(nSeqsBreaks[process], nSeqsBreaks[process+1], child_otuIndex[process], child_seqIndex[process], child_singleTau[process]);
1626                                         totals[process] = child_otuIndex[process].size();
1627
1628                                         exit(0);
1629                                 }else { 
1630                                         m->mothurOut("[ERROR]: unable to spawn the necessary processes. Error code: " + toString(pid)); m->mothurOutEndLine(); 
1631                                         perror(" : ");
1632                                         for (int i=0;i<processIDs.size();i++) {  int temp = processIDs[i]; kill (temp, SIGINT); }
1633                                         exit(0);
1634                                 }
1635                         }
1636                         
1637                         //parent does its part
1638                         calcNewDistancesParent(nSeqsBreaks[0], nSeqsBreaks[1]);
1639                         int total = seqIndex.size();
1640
1641                         //force parent to wait until all the processes are done
1642                         for (int i=0;i<processIDs.size();i++) { 
1643                                 int temp = processIDs[i];
1644                                 wait(&temp);
1645                         }
1646
1647                         for(int i=1;i<processors;i++){
1648                                 int oldTotal = total;
1649                                 total += totals[i];
1650
1651                                 singleTau.resize(total, 0);
1652                                 seqIndex.resize(total, 0);
1653                                 seqNumber.resize(total, 0);
1654                                 
1655                                 int childIndex = 0;
1656                                 
1657                                 for(int j=oldTotal;j<total;j++){
1658                                         int otuI = child_otuIndex[i][childIndex];
1659                                         int seqI = child_seqIndex[i][childIndex];
1660
1661                                         singleTau[j] = child_singleTau[i][childIndex];
1662                                         aaP[otuI][nSeqsPerOTU[otuI]] = j;
1663                                         aaI[otuI][nSeqsPerOTU[otuI]] = seqI;
1664                                         nSeqsPerOTU[otuI]++;
1665
1666                                         childIndex++;
1667                                 }
1668                         }
1669                 }
1670 #else
1671                 calcNewDistancesParent(0, numSeqs);             
1672 #endif          
1673         }
1674         catch(exception& e) {
1675                 m->errorOut(e, "ShhherCommand", "calcNewDistances");
1676                 exit(1);        
1677         }               
1678 }
1679
1680 /**************************************************************************************************/
1681 #ifdef USE_MPI
1682 void ShhherCommand::calcNewDistancesChildMPI(int startSeq, int stopSeq, vector<int>& otuIndex){
1683         
1684         try{
1685                 vector<double> newTau(numOTUs,0);
1686                 vector<double> norms(numSeqs, 0);
1687                 otuIndex.clear();
1688                 seqIndex.clear();
1689                 singleTau.clear();
1690                 
1691                 
1692                 
1693                 for(int i=startSeq;i<stopSeq;i++){
1694                         double offset = 1e8;
1695                         int indexOffset = i * numOTUs;
1696                         
1697                         for(int j=0;j<numOTUs;j++){
1698                                 
1699                                 if(weight[j] > MIN_WEIGHT && change[j] == 1){
1700                                         dist[indexOffset + j] = getDistToCentroid(centroids[j], i, lengths[i]);
1701                                 }
1702                                 if(weight[j] > MIN_WEIGHT && dist[indexOffset + j] < offset){
1703                                         offset = dist[indexOffset + j];
1704                                 }
1705                         }
1706                         
1707                         for(int j=0;j<numOTUs;j++){
1708                                 if(weight[j] > MIN_WEIGHT){
1709                                         newTau[j] = exp(sigma * (-dist[indexOffset + j] + offset)) * weight[j];
1710                                         norms[i] += newTau[j];
1711                                 }
1712                                 else{
1713                                         newTau[j] = 0.0;
1714                                 }
1715                         }
1716                         
1717                         for(int j=0;j<numOTUs;j++){
1718
1719                                 newTau[j] /= norms[i];
1720                                 
1721                                 if(newTau[j] > MIN_TAU){
1722                                         otuIndex.push_back(j);
1723                                         seqIndex.push_back(i);
1724                                         singleTau.push_back(newTau[j]);
1725                                 }
1726                         }
1727                         
1728                 }
1729         }
1730         catch(exception& e) {
1731                 m->errorOut(e, "ShhherCommand", "calcNewDistancesChildMPI");
1732                 exit(1);        
1733         }               
1734 }
1735 #endif
1736 /**************************************************************************************************/
1737
1738 void ShhherCommand::calcNewDistancesChild(int startSeq, int stopSeq, vector<int>& child_otuIndex, vector<int>& child_seqIndex, vector<double>& child_singleTau){
1739         
1740         try{
1741                 vector<double> newTau(numOTUs,0);
1742                 vector<double> norms(numSeqs, 0);
1743                 child_otuIndex.resize(0);
1744                 child_seqIndex.resize(0);
1745                 child_singleTau.resize(0);
1746                 
1747                 for(int i=startSeq;i<stopSeq;i++){
1748                         double offset = 1e8;
1749                         int indexOffset = i * numOTUs;
1750                         
1751                         
1752                         for(int j=0;j<numOTUs;j++){
1753                                 if(weight[j] > MIN_WEIGHT && change[j] == 1){
1754                                         dist[indexOffset + j] = getDistToCentroid(centroids[j], i, lengths[i]);
1755                                 }
1756                                 
1757                                 if(weight[j] > MIN_WEIGHT && dist[indexOffset + j] < offset){
1758                                         offset = dist[indexOffset + j];
1759                                 }
1760                         }
1761                         
1762                         for(int j=0;j<numOTUs;j++){
1763                                 if(weight[j] > MIN_WEIGHT){
1764                                         newTau[j] = exp(sigma * (-dist[indexOffset + j] + offset)) * weight[j];
1765                                         norms[i] += newTau[j];
1766                                 }
1767                                 else{
1768                                         newTau[j] = 0.0;
1769                                 }
1770                         }
1771                         
1772                         for(int j=0;j<numOTUs;j++){
1773                                 newTau[j] /= norms[i];
1774                                 
1775                                 if(newTau[j] > MIN_TAU){
1776                                         child_otuIndex.push_back(j);
1777                                         child_seqIndex.push_back(i);
1778                                         child_singleTau.push_back(newTau[j]);
1779                                 }
1780                         }
1781                 }
1782         }
1783         catch(exception& e) {
1784                 m->errorOut(e, "ShhherCommand", "calcNewDistancesChild");
1785                 exit(1);        
1786         }               
1787 }
1788
1789 /**************************************************************************************************/
1790
1791 void ShhherCommand::calcNewDistancesParent(int startSeq, int stopSeq){
1792         
1793         try{
1794                 
1795                 int total = 0;
1796                 vector<double> newTau(numOTUs,0);
1797                 vector<double> norms(numSeqs, 0);
1798                 nSeqsPerOTU.assign(numOTUs, 0);
1799                 
1800                 for(int i=startSeq;i<stopSeq;i++){
1801                         int indexOffset = i * numOTUs;
1802                         
1803                         double offset = 1e8;
1804                         
1805                         for(int j=0;j<numOTUs;j++){
1806                                 if(weight[j] > MIN_WEIGHT && change[j] == 1){
1807                                         dist[indexOffset + j] = getDistToCentroid(centroids[j], i, lengths[i]);
1808                                 }
1809                                 
1810                                 if(weight[j] > MIN_WEIGHT && dist[indexOffset + j] < offset){
1811                                         offset = dist[indexOffset + j];
1812                                 }
1813                         }
1814                         
1815                         for(int j=0;j<numOTUs;j++){
1816                                 if(weight[j] > MIN_WEIGHT){
1817                                         newTau[j] = exp(sigma * (-dist[indexOffset + j] + offset)) * weight[j];
1818                                         norms[i] += newTau[j];
1819                                 }
1820                                 else{
1821                                         newTau[j] = 0.0;
1822                                 }
1823                         }
1824                         
1825                         for(int j=0;j<numOTUs;j++){
1826                                 newTau[j] /= norms[i];
1827                         }
1828                         
1829                         for(int j=0;j<numOTUs;j++){
1830                                 if(newTau[j] > MIN_TAU){
1831                                         
1832                                         int oldTotal = total;
1833                                         
1834                                         total++;
1835                                         
1836                                         singleTau.resize(total, 0);
1837                                         seqNumber.resize(total, 0);
1838                                         seqIndex.resize(total, 0);
1839                                         
1840                                         singleTau[oldTotal] = newTau[j];
1841                                         
1842                                         aaP[j][nSeqsPerOTU[j]] = oldTotal;
1843                                         aaI[j][nSeqsPerOTU[j]] = i;
1844                                         nSeqsPerOTU[j]++;
1845                                 }
1846                         }
1847                 }
1848         }
1849         catch(exception& e) {
1850                 m->errorOut(e, "ShhherCommand", "calcNewDistancesParent");
1851                 exit(1);        
1852         }               
1853 }
1854
1855 /**************************************************************************************************/
1856
1857 void ShhherCommand::setOTUs(){
1858         
1859         try {
1860                 vector<double> bigTauMatrix(numOTUs * numSeqs, 0.0000);
1861                 
1862                 for(int i=0;i<numOTUs;i++){
1863                         for(int j=0;j<nSeqsPerOTU[i];j++){
1864                                 int index = cumNumSeqs[i] + j;
1865                                 double tauValue = singleTau[seqNumber[index]];
1866                                 int sIndex = seqIndex[index];
1867                                 bigTauMatrix[sIndex * numOTUs + i] = tauValue;                          
1868                         }
1869                 }
1870                 
1871                 for(int i=0;i<numSeqs;i++){
1872                         double maxTau = -1.0000;
1873                         int maxOTU = -1;
1874                         for(int j=0;j<numOTUs;j++){
1875                                 if(bigTauMatrix[i * numOTUs + j] > maxTau){
1876                                         maxTau = bigTauMatrix[i * numOTUs + j];
1877                                         maxOTU = j;
1878                                 }
1879                         }
1880                         
1881                         otuData[i] = maxOTU;
1882                 }
1883                 
1884                 nSeqsPerOTU.assign(numOTUs, 0);         
1885                 
1886                 for(int i=0;i<numSeqs;i++){
1887                         int index = otuData[i];
1888                         
1889                         singleTau[i] = 1.0000;
1890                         dist[i] = 0.0000;
1891                         
1892                         aaP[index][nSeqsPerOTU[index]] = i;
1893                         aaI[index][nSeqsPerOTU[index]] = i;
1894                         
1895                         nSeqsPerOTU[index]++;
1896                 }
1897                 fill(); 
1898         }
1899         catch(exception& e) {
1900                 m->errorOut(e, "ShhherCommand", "calcNewDistances");
1901                 exit(1);        
1902         }               
1903 }
1904
1905 /**************************************************************************************************/
1906
1907 void ShhherCommand::writeQualities(vector<int> otuCounts){
1908         
1909         try {
1910                 string qualityFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.qual";
1911
1912                 ofstream qualityFile;
1913                 m->openOutputFile(qualityFileName, qualityFile);
1914
1915                 qualityFile.setf(ios::fixed, ios::floatfield);
1916                 qualityFile.setf(ios::showpoint);
1917                 qualityFile << setprecision(6);
1918                 
1919                 vector<vector<int> > qualities(numOTUs);
1920                 vector<double> pr(HOMOPS, 0);
1921                 
1922                 
1923                 for(int i=0;i<numOTUs;i++){
1924                         int index = 0;
1925                         int base = 0;
1926                         
1927                         if(nSeqsPerOTU[i] > 0){
1928                                 qualities[i].assign(1024, -1);
1929                                 
1930                                 while(index < numFlowCells){
1931                                         double maxPrValue = 1e8;
1932                                         short maxPrIndex = -1;
1933                                         double count = 0.0000;
1934                                         
1935                                         pr.assign(HOMOPS, 0);
1936                                         
1937                                         for(int j=0;j<nSeqsPerOTU[i];j++){
1938                                                 int lIndex = cumNumSeqs[i] + j;
1939                                                 double tauValue = singleTau[seqNumber[lIndex]];
1940                                                 int sequenceIndex = aaI[i][j];
1941                                                 short intensity = flowDataIntI[sequenceIndex * numFlowCells + index];
1942                                                 
1943                                                 count += tauValue;
1944                                                 
1945                                                 for(int s=0;s<HOMOPS;s++){
1946                                                         pr[s] += tauValue * singleLookUp[s * NUMBINS + intensity];
1947                                                 }
1948                                         }
1949                                         
1950                                         maxPrIndex = uniqueFlowgrams[centroids[i] * numFlowCells + index];
1951                                         maxPrValue = pr[maxPrIndex];
1952                                         
1953                                         if(count > MIN_COUNT){
1954                                                 double U = 0.0000;
1955                                                 double norm = 0.0000;
1956                                                 
1957                                                 for(int s=0;s<HOMOPS;s++){
1958                                                         norm += exp(-(pr[s] - maxPrValue));
1959                                                 }
1960                                                 
1961                                                 for(int s=1;s<=maxPrIndex;s++){
1962                                                         int value = 0;
1963                                                         double temp = 0.0000;
1964                                                         
1965                                                         U += exp(-(pr[s-1]-maxPrValue))/norm;
1966                                                         
1967                                                         if(U>0.00){
1968                                                                 temp = log10(U);
1969                                                         }
1970                                                         else{
1971                                                                 temp = -10.1;
1972                                                         }
1973                                                         temp = floor(-10 * temp);
1974                                                         value = (int)floor(temp);
1975                                                         if(value > 100){        value = 100;    }
1976                                                         
1977                                                         qualities[i][base] = (int)value;
1978                                                         base++;
1979                                                 }
1980                                         }
1981                                         
1982                                         index++;
1983                                 }
1984                         }
1985                         
1986                         
1987                         if(otuCounts[i] > 0){
1988                                 qualityFile << '>' << seqNameVector[mapUniqueToSeq[i]] << endl;
1989                                 
1990                                 int j=4;        //need to get past the first four bases
1991                                 while(qualities[i][j] != -1){
1992                                         qualityFile << qualities[i][j] << ' ';
1993                                         j++;
1994                                 }
1995                                 qualityFile << endl;
1996                         }
1997                 }
1998                 qualityFile.close();
1999                 
2000         }
2001         catch(exception& e) {
2002                 m->errorOut(e, "ShhherCommand", "writeQualities");
2003                 exit(1);        
2004         }               
2005 }
2006
2007 /**************************************************************************************************/
2008
2009 void ShhherCommand::writeSequences(vector<int> otuCounts){
2010         try {
2011                 string bases = "TACG";
2012                 
2013                 string fastaFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.fasta";
2014                 ofstream fastaFile;
2015                 m->openOutputFile(fastaFileName, fastaFile);
2016                 
2017                 vector<string> names(numOTUs, "");
2018                 
2019                 for(int i=0;i<numOTUs;i++){
2020                         int index = centroids[i];
2021                         
2022                         if(otuCounts[i] > 0){
2023                                 fastaFile << '>' << seqNameVector[aaI[i][0]] << endl;
2024                                 
2025                                 for(int j=8;j<numFlowCells;j++){
2026                                         
2027                                         char base = bases[j % 4];
2028                                         for(int k=0;k<uniqueFlowgrams[index * numFlowCells + j];k++){
2029                                                 fastaFile << base;
2030                                         }
2031                                 }
2032                                 fastaFile << endl;
2033                         }
2034                 }
2035                 fastaFile.close();
2036                 
2037                 if(compositeFASTAFileName != ""){
2038                         m->appendFiles(fastaFileName, compositeFASTAFileName);
2039                 }
2040         }
2041         catch(exception& e) {
2042                 m->errorOut(e, "ShhherCommand", "writeSequences");
2043                 exit(1);        
2044         }               
2045 }
2046
2047 /**************************************************************************************************/
2048
2049 void ShhherCommand::writeNames(vector<int> otuCounts){
2050         try {
2051                 string nameFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.final.names";
2052                 ofstream nameFile;
2053                 m->openOutputFile(nameFileName, nameFile);
2054                 
2055                 for(int i=0;i<numOTUs;i++){
2056                         if(otuCounts[i] > 0){
2057                                 nameFile << seqNameVector[aaI[i][0]] << '\t' << seqNameVector[aaI[i][0]];
2058                                 
2059                                 for(int j=1;j<nSeqsPerOTU[i];j++){
2060                                         nameFile << ',' << seqNameVector[aaI[i][j]];
2061                                 }
2062                                 
2063                                 nameFile << endl;
2064                         }
2065                 }
2066                 nameFile.close();
2067         }
2068         catch(exception& e) {
2069                 m->errorOut(e, "ShhherCommand", "writeNames");
2070                 exit(1);        
2071         }               
2072 }
2073
2074 /**************************************************************************************************/
2075
2076 void ShhherCommand::writeGroups(){
2077         try {
2078                 string fileRoot = flowFileName.substr(0,flowFileName.find_last_of('.'));
2079                 string groupFileName = fileRoot + ".pn.groups";
2080                 ofstream groupFile;
2081                 m->openOutputFile(groupFileName, groupFile);
2082                 
2083                 for(int i=0;i<numSeqs;i++){
2084                         groupFile << seqNameVector[i] << '\t' << fileRoot << endl;
2085                 }
2086                 groupFile.close();
2087         }
2088         catch(exception& e) {
2089                 m->errorOut(e, "ShhherCommand", "writeGroups");
2090                 exit(1);        
2091         }               
2092 }
2093
2094 /**************************************************************************************************/
2095
2096 void ShhherCommand::writeClusters(vector<int> otuCounts){
2097         try {
2098                 string otuCountsFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.counts";
2099                 ofstream otuCountsFile;
2100                 m->openOutputFile(otuCountsFileName, otuCountsFile);
2101                 
2102                 string bases = "TACG";
2103                 
2104                 for(int i=0;i<numOTUs;i++){
2105                         //output the translated version of the centroid sequence for the otu
2106                         if(otuCounts[i] > 0){
2107                                 int index = centroids[i];
2108                                 
2109                                 otuCountsFile << "ideal\t";
2110                                 for(int j=8;j<numFlowCells;j++){
2111                                         char base = bases[j % 4];
2112                                         for(int s=0;s<uniqueFlowgrams[index * numFlowCells + j];s++){
2113                                                 otuCountsFile << base;
2114                                         }
2115                                 }
2116                                 otuCountsFile << endl;
2117                                 
2118                                 for(int j=0;j<nSeqsPerOTU[i];j++){
2119                                         int sequence = aaI[i][j];
2120                                         otuCountsFile << seqNameVector[sequence] << '\t';
2121                                         
2122                                         for(int k=8;k<lengths[sequence];k++){
2123                                                 char base = bases[k % 4];
2124                                                 int freq = int(0.01 * (double)flowDataIntI[sequence * numFlowCells + k] + 0.5);
2125                                                 
2126                                                 for(int s=0;s<freq;s++){
2127                                                         otuCountsFile << base;
2128                                                 }
2129                                         }
2130                                         otuCountsFile << endl;
2131                                 }
2132                                 otuCountsFile << endl;
2133                         }
2134                 }
2135                 otuCountsFile.close();
2136         }
2137         catch(exception& e) {
2138                 m->errorOut(e, "ShhherCommand", "writeClusters");
2139                 exit(1);        
2140         }               
2141 }
2142
2143 //**********************************************************************************************************************