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