]> git.donarmstrong.com Git - rsem.git/blob - EM.cpp
RSEM Source Codes
[rsem.git] / EM.cpp
1 #include<ctime>
2 #include<cmath>
3 #include<cstdio>
4 #include<cstdlib>
5 #include<cstring>
6 #include<cassert>
7 #include<string>
8 #include<vector>
9 #include<algorithm>
10 #include<pthread.h>
11
12 #include "utils.h"
13
14 #include "Read.h"
15 #include "SingleRead.h"
16 #include "SingleReadQ.h"
17 #include "PairedEndRead.h"
18 #include "PairedEndReadQ.h"
19
20 #include "SingleHit.h"
21 #include "PairedEndHit.h"
22
23 #include "Model.h"
24 #include "SingleModel.h"
25 #include "SingleQModel.h"
26 #include "PairedEndModel.h"
27 #include "PairedEndQModel.h"
28
29 #include "Transcript.h"
30 #include "Transcripts.h"
31
32 #include "Refs.h"
33 #include "GroupInfo.h"
34 #include "HitContainer.h"
35 #include "ReadIndex.h"
36 #include "ReadReader.h"
37
38 #include "ModelParams.h"
39
40 #include "HitWrapper.h"
41 #include "BamWriter.h"
42
43 using namespace std;
44
45 const double STOP_CRITERIA = 0.001;
46 const int MAX_ROUND = 10000;
47 const int MIN_ROUND = 20;
48
49 struct Params {
50         void *model;
51         void *reader, *hitv, *ncpv, *mhp, *countv;
52 };
53
54 int read_type;
55 int m, M; // m genes, M isoforms
56 int N0, N1, N2, N_tot;
57 int nThreads;
58
59
60 bool genBamF; // If user wants to generate bam file, true; otherwise, false.
61 bool updateModel, calcExpectedWeights;
62 bool genGibbsOut; // generate file for Gibbs sampler
63
64 char refName[STRLEN], imdName[STRLEN], outName[STRLEN];
65 char refF[STRLEN], groupF[STRLEN], cntF[STRLEN], tiF[STRLEN];
66 char mparamsF[STRLEN], bmparamsF[STRLEN];
67
68 char inpSamType;
69 char *pt_fn_list, *pt_chr_list;
70 char inpSamF[STRLEN], outBamF[STRLEN], fn_list[STRLEN], chr_list[STRLEN];
71
72 char out_for_gibbs_F[STRLEN];
73
74 vector<double> theta, eel; // eel : expected effective length
75
76 double *probv, **countvs;
77
78 Refs refs;
79 GroupInfo gi;
80 Transcripts transcripts;
81
82 ModelParams mparams;
83
84 template<class ReadType, class HitType, class ModelType>
85 void init(ReadReader<ReadType> **&readers, HitContainer<HitType> **&hitvs, double **&ncpvs, ModelType **&mhps) {
86         int nReads, nHits, rt;
87         int nrLeft, nhT, curnr; // nrLeft : number of reads left, nhT : hit threshold per thread, curnr: current number of reads
88         char datF[STRLEN];
89
90         int s;
91         char readFs[2][STRLEN];
92         ReadIndex *indices[2];
93         ifstream fin;
94
95         readers = new ReadReader<ReadType>*[nThreads];
96         genReadFileNames(imdName, 1, read_type, s, readFs);
97         for (int i = 0; i < s; i++) {
98                 indices[i] = new ReadIndex(readFs[i]);
99         }
100         for (int i = 0; i < nThreads; i++) {
101                 readers[i] = new ReadReader<ReadType>(s, readFs);
102                 readers[i]->setIndices(indices);
103         }
104
105         hitvs = new HitContainer<HitType>*[nThreads];
106         for (int i = 0; i < nThreads; i++) {
107                 hitvs[i] = new HitContainer<HitType>();
108         }
109
110         sprintf(datF, "%s.dat", imdName);
111         fin.open(datF);
112         if (!fin.is_open()) { fprintf(stderr, "Cannot open %s! It may not exist.\n", datF); exit(-1); }
113         fin>>nReads>>nHits>>rt;
114         if (nReads != N1) { fprintf(stderr, "Number of alignable reads does not match!\n"); exit(-1); }
115         //assert(nReads == N1);
116         if (rt != read_type) { fprintf(stderr, "Data file (.dat) does not have the right read type!\n"); exit(-1); }
117         //assert(rt == read_type);
118
119         //A just so so strategy for paralleling
120         nhT = nHits / nThreads;
121         nrLeft = N1;
122         curnr = 0;
123
124         ncpvs = new double*[nThreads];
125         for (int i = 0; i < nThreads; i++) {
126                 int ntLeft = nThreads - i - 1; // # of threads left
127                 if (!readers[i]->locate(curnr)) { fprintf(stderr, "Read indices files do not match!\n"); exit(-1); }
128                 //assert(readers[i]->locate(curnr));
129
130                 while (nrLeft > ntLeft && hitvs[i]->getNHits() < nhT) {
131                         if (!hitvs[i]->read(fin)) { fprintf(stderr, "Cannot read alignments from .dat file!\n"); exit(-1); }
132                         //assert(hitvs[i]->read(fin));
133                         --nrLeft;
134                         if (verbose && nrLeft % 1000000 == 0) { printf("DAT %d reads left!\n", nrLeft); }
135                 }
136                 ncpvs[i] = new double[hitvs[i]->getN()];
137                 memset(ncpvs[i], 0, sizeof(double) * hitvs[i]->getN());
138                 curnr += hitvs[i]->getN();
139
140                 if (verbose) { printf("Thread %d : N = %d, NHit = %d\n", i, hitvs[i]->getN(), hitvs[i]->getNHits()); }
141         }
142
143         fin.close();
144
145         mhps = new ModelType*[nThreads];
146         for (int i = 0; i < nThreads; i++) {
147                 mhps[i] = new ModelType(mparams, false); // just model helper
148         }
149
150         probv = new double[M + 1];
151         countvs = new double*[nThreads];
152         for (int i = 0; i < nThreads; i++) {
153                 countvs[i] = new double[M + 1];
154         }
155
156
157         if (verbose) { printf("EM_init finished!\n"); }
158 }
159
160 template<class ReadType, class HitType, class ModelType>
161 void* E_STEP(void* arg) {
162         Params *params = (Params*)arg;
163         ModelType *model = (ModelType*)(params->model);
164         ReadReader<ReadType> *reader = (ReadReader<ReadType>*)(params->reader);
165         HitContainer<HitType> *hitv = (HitContainer<HitType>*)(params->hitv);
166         double *ncpv = (double*)(params->ncpv);
167         ModelType *mhp = (ModelType*)(params->mhp);
168         double *countv = (double*)(params->countv);
169
170         bool needCalcConPrb = model->getNeedCalcConPrb();
171
172         ReadType read;
173
174         int N = hitv->getN();
175         double sum;
176         vector<double> fracs; //to remove this, do calculation twice
177         int fr, to, id;
178
179         if (needCalcConPrb || updateModel) { reader->reset(); }
180         if (updateModel) { mhp->init(); }
181
182         memset(countv, 0, sizeof(double) * (M + 1));
183         for (int i = 0; i < N; i++) {
184                 if (needCalcConPrb || updateModel) {
185                         if (!reader->next(read)) {
186                                 fprintf(stderr, "Can not load a read!\n");
187                                 exit(-1);
188                         }
189                         //assert(reader->next(read));
190                 }
191                 fr = hitv->getSAt(i);
192                 to = hitv->getSAt(i + 1);
193                 fracs.resize(to - fr + 1);
194
195                 sum = 0.0;
196
197                 if (needCalcConPrb) { ncpv[i] = model->getNoiseConPrb(read); }
198                 fracs[0] = probv[0] * ncpv[i];
199                 if (fracs[0] < EPSILON) fracs[0] = 0.0;
200                 sum += fracs[0];
201                 for (int j = fr; j < to; j++) {
202                         HitType &hit = hitv->getHitAt(j);
203                         if (needCalcConPrb) { hit.setConPrb(model->getConPrb(read, hit)); }
204                         id = j - fr + 1;
205                         fracs[id] = probv[hit.getSid()] * hit.getConPrb();
206                         if (fracs[id] < EPSILON) fracs[id] = 0.0;
207                         sum += fracs[id];
208                 }
209
210                 if (sum >= EPSILON) {
211                         fracs[0] /= sum;
212                         countv[0] += fracs[0];
213                         if (updateModel) { mhp->updateNoise(read, fracs[0]); }
214                         if (calcExpectedWeights) { ncpv[i] = fracs[0]; }
215                         for (int j = fr; j < to; j++) {
216                                 HitType &hit = hitv->getHitAt(j);
217                                 id = j - fr + 1;
218                                 fracs[id] /= sum;
219                                 countv[hit.getSid()] += fracs[id];
220                                 if (updateModel) { mhp->update(read, hit, fracs[id]); }
221                                 if (calcExpectedWeights) { hit.setConPrb(fracs[id]); }
222                         }                       
223                 }
224                 else if (calcExpectedWeights) {
225                         ncpv[i] = 0.0;
226                         for (int j = fr; j < to; j++) {
227                                 HitType &hit = hitv->getHitAt(j);
228                                 hit.setConPrb(0.0);
229                         }
230                 }
231         }
232
233         return NULL;
234 }
235
236 template<class ReadType, class HitType, class ModelType>
237 void* calcConProbs(void* arg) {
238         Params *params = (Params*)arg;
239         ModelType *model = (ModelType*)(params->model);
240         ReadReader<ReadType> *reader = (ReadReader<ReadType>*)(params->reader);
241         HitContainer<HitType> *hitv = (HitContainer<HitType>*)(params->hitv);
242         double *ncpv = (double*)(params->ncpv);
243
244         ReadType read;
245         int N = hitv->getN();
246         int fr, to;
247
248         assert(model->getNeedCalcConPrb());
249         reader->reset();
250
251         for (int i = 0; i < N; i++) {
252                 if (!reader->next(read)) {
253                         fprintf(stderr, "Can not load a read!\n");
254                         exit(-1);
255                 }
256                 fr = hitv->getSAt(i);
257                 to = hitv->getSAt(i + 1);
258
259                 ncpv[i] = model->getNoiseConPrb(read);
260                 for (int j = fr; j < to; j++) {
261                         HitType &hit = hitv->getHitAt(j);
262                         hit.setConPrb(model->getConPrb(read, hit));
263                 }
264         }
265
266         return NULL;
267 }
268
269 template<class ModelType>
270 void calcExpectedEffectiveLengths(ModelType& model) {
271   int lb, ub, span;
272   double *pdf = NULL, *cdf = NULL, *clen = NULL; // clen[i] = sigma_{j=1}^{i}pdf[i]*(lb+i)
273   
274   model.getGLD().copyTo(pdf, cdf, lb, ub, span);
275   clen = new double[span + 1];
276   clen[0] = 0.0;
277   for (int i = 1; i <= span; i++) {
278     clen[i] = clen[i - 1] + pdf[i] * (lb + i);
279   }
280
281   eel.clear();
282   eel.resize(M + 1, 0.0);
283   for (int i = 1; i <= M; i++) {
284     int totLen = refs.getRef(i).getTotLen();
285     int fullLen = refs.getRef(i).getFullLen();
286     int pos1 = max(min(totLen - fullLen + 1, ub) - lb, 0);
287     int pos2 = max(min(totLen, ub) - lb, 0);
288
289     if (pos2 == 0) { eel[i] = 0.0; continue; }
290     
291     eel[i] = fullLen * cdf[pos1] + ((cdf[pos2] - cdf[pos1]) * (totLen + 1) - (clen[pos2] - clen[pos1]));
292     assert(eel[i] >= 0);
293     if (eel[i] < MINEEL) { eel[i] = 0.0; }
294   }
295   
296   delete[] pdf;
297   delete[] cdf;
298   delete[] clen;
299 }
300
301 template<class ModelType>
302 void writeResults(ModelType& model, double* counts) {
303         double denom;
304         char modelF[STRLEN], thetaF[STRLEN];
305         char outF[STRLEN];
306         FILE *fo;
307
308         sprintf(modelF, "%s.model", outName);
309         model.write(modelF);
310
311         sprintf(thetaF, "%s.theta", outName);
312         fo = fopen(thetaF, "w");
313         fprintf(fo, "%d\n", M + 1);
314         for (int i = 0; i < M; i++) fprintf(fo, "%.15g ", theta[i]);
315         fprintf(fo, "%.15g\n", theta[M]);
316         fclose(fo);
317
318
319         //calculate normalized read fraction
320         double *nrf = new double[M + 1];
321         memset(nrf, 0, sizeof(double) * (M + 1));
322         denom = 1.0 - theta[0];
323         if (denom <= 0) { fprintf(stderr, "No alignable reads?!\n"); exit(-1); }
324         for (int i = 1; i <= M; i++) nrf[i] = theta[i] / denom;
325
326         //calculate tau values
327         double *tau = new double[M + 1];
328         memset(tau, 0, sizeof(double) * (M + 1));
329
330         denom = 0.0;
331         for (int i = 1; i <= M; i++) 
332           if (eel[i] > EPSILON) {
333             tau[i] = theta[i] / eel[i];
334             denom += tau[i];
335           }   
336         if (denom <= 0) { fprintf(stderr, "No alignable reads?!\n"); exit(-1); }
337         //assert(denom > 0);
338         for (int i = 1; i <= M; i++) {
339                 tau[i] /= denom;
340         }
341
342         //isoform level results
343         sprintf(outF, "%s.iso_res", imdName);
344         fo = fopen(outF, "w");
345         for (int i = 1; i <= M; i++) {
346                 const Transcript& transcript = transcripts.getTranscriptAt(i);
347                 fprintf(fo, "%s%c", transcript.getTranscriptID().c_str(), (i < M ? '\t' : '\n'));
348         }
349         for (int i = 1; i <= M; i++)
350                 fprintf(fo, "%.2f%c", counts[i], (i < M ? '\t' : '\n'));
351         for (int i = 1; i <= M; i++)
352                 fprintf(fo, "%.15g%c", nrf[i], (i < M ? '\t' : '\n'));
353         for (int i = 1; i <= M; i++)
354                 fprintf(fo, "%.15g%c", tau[i], (i < M ? '\t' : '\n'));
355         for (int i = 1; i <= M; i++) {
356                 const Transcript& transcript = transcripts.getTranscriptAt(i);
357                 fprintf(fo, "%s%c", transcript.getLeft().c_str(), (i < M ? '\t' : '\n'));
358         }
359         fclose(fo);
360
361         //gene level results
362         sprintf(outF, "%s.gene_res", imdName);
363         fo = fopen(outF, "w");
364         for (int i = 0; i < m; i++) {
365                 const string& gene_id = transcripts.getTranscriptAt(gi.spAt(i)).getGeneID();
366                 fprintf(fo, "%s%c", gene_id.c_str(), (i < m - 1 ? '\t' : '\n'));
367         }
368         for (int i = 0; i < m; i++) {
369                 double sumC = 0.0; // sum of counts
370                 int b = gi.spAt(i), e = gi.spAt(i + 1);
371                 for (int j = b; j < e; j++) sumC += counts[j];
372                 fprintf(fo, "%.2f%c", sumC, (i < m - 1 ? '\t' : '\n'));
373         }
374         for (int i = 0; i < m; i++) {
375                 double sumN = 0.0; // sum of normalized read fraction
376                 int b = gi.spAt(i), e = gi.spAt(i + 1);
377                 for (int j = b; j < e; j++) sumN += nrf[j];
378                 fprintf(fo, "%.15g%c", sumN, (i < m - 1 ? '\t' : '\n'));
379         }
380         for (int i = 0; i < m; i++) {
381                 double sumT = 0.0; // sum of tau values
382                 int b = gi.spAt(i), e = gi.spAt(i + 1);
383                 for (int j = b; j < e; j++) sumT += tau[j];
384                 fprintf(fo, "%.15g%c", sumT, (i < m - 1 ? '\t' : '\n'));
385         }
386         for (int i = 0; i < m; i++) {
387                 int b = gi.spAt(i), e = gi.spAt(i + 1);
388                 for (int j = b; j < e; j++) {
389                         fprintf(fo, "%s%c", transcripts.getTranscriptAt(j).getTranscriptID().c_str(), (j < e - 1 ? ',' : (i < m - 1 ? '\t' :'\n')));
390                 }
391         }
392         fclose(fo);
393
394         delete[] nrf;
395         delete[] tau;
396
397         if (verbose) { printf("Expression Results are written!\n"); }
398 }
399
400 template<class ReadType, class HitType, class ModelType>
401 void release(ReadReader<ReadType> **readers, HitContainer<HitType> **hitvs, double **ncpvs, ModelType **mhps) {
402         delete[] probv;
403         for (int i = 0; i < nThreads; i++) {
404                 delete[] countvs[i];
405         }
406         delete[] countvs;
407
408         for (int i = 0; i < nThreads; i++) {
409                 delete readers[i];
410                 delete hitvs[i];
411                 delete[] ncpvs[i];
412                 delete mhps[i];
413         }
414         delete[] readers;
415         delete[] hitvs;
416         delete[] ncpvs;
417         delete[] mhps;
418 }
419
420 inline bool doesUpdateModel(int ROUND) {
421         //return false; // never update, for debugging only
422         return ROUND <= 20 || ROUND % 100 == 0;
423 }
424
425 //Including initialize, algorithm and results saving
426 template<class ReadType, class HitType, class ModelType>
427 void EM() {
428         int ROUND;
429         double sum;
430
431         double bChange = 0.0, change = 0.0; // bChange : biggest change
432         int totNum = 0;
433
434         ModelType model(mparams); //master model
435         ReadReader<ReadType> **readers;
436         HitContainer<HitType> **hitvs;
437         double **ncpvs;
438         ModelType **mhps; //model helpers
439
440         Params fparams[nThreads];
441         pthread_t threads[nThreads];
442         pthread_attr_t attr;
443         void *status;
444         int rc;
445
446         //initialize boolean variables
447         updateModel = calcExpectedWeights = false;
448
449         theta.clear();
450         theta.resize(M + 1, 0.0);
451         init<ReadType, HitType, ModelType>(readers, hitvs, ncpvs, mhps);
452
453         //set initial parameters
454         assert(N_tot > N2);
455         theta[0] = max(N0 * 1.0 / (N_tot - N2), 1e-8);
456         double val = (1.0 - theta[0]) / M;
457         for (int i = 1; i <= M; i++) theta[i] = val;
458
459         model.estimateFromReads(imdName);
460
461         for (int i = 0; i < nThreads; i++) {
462                 fparams[i].model = (void*)(&model);
463
464                 fparams[i].reader = (void*)readers[i];
465                 fparams[i].hitv = (void*)hitvs[i];
466                 fparams[i].ncpv = (void*)ncpvs[i];
467                 fparams[i].mhp = (void*)mhps[i];
468                 fparams[i].countv = (void*)countvs[i];
469         }
470
471         /* set thread attribute to be joinable */
472         pthread_attr_init(&attr);
473         pthread_attr_setdetachstate(&attr, PTHREAD_CREATE_JOINABLE);
474
475
476         ROUND = 0;
477         do {
478                 ++ROUND;
479
480                 updateModel = doesUpdateModel(ROUND);
481
482                 for (int i = 0; i <= M; i++) probv[i] = theta[i];
483
484                 //E step
485                 for (int i = 0; i < nThreads; i++) {
486                         rc = pthread_create(&threads[i], &attr, E_STEP<ReadType, HitType, ModelType>, (void*)(&fparams[i]));
487                         if (rc != 0) { fprintf(stderr, "Cannot create thread %d at ROUND %d! (numbered from 0)\n", i, ROUND); exit(-1); }
488                         //assert(rc == 0);
489                 }
490
491                 for (int i = 0; i < nThreads; i++) {
492                         rc = pthread_join(threads[i], &status);
493                         if (rc != 0) { fprintf(stderr, "Cannot join thread %d at ROUND %d! (numbered from 0)\n", i, ROUND); exit(-1); }
494                         //assert(rc == 0);
495                 }
496
497                 model.setNeedCalcConPrb(false);
498
499                 for (int i = 1; i < nThreads; i++) {
500                         for (int j = 0; j <= M; j++) {
501                                 countvs[0][j] += countvs[i][j];
502                         }
503                 }
504
505                 //add N0 noise reads
506                 countvs[0][0] += N0;
507
508                 //M step;
509                 sum = 0.0;
510                 for (int i = 0; i <= M; i++) sum += countvs[0][i];
511                 assert(sum >= EPSILON);
512                 for (int i = 0; i <= M; i++) theta[i] = countvs[0][i] / sum;
513
514                 if (updateModel) {
515                         model.init();
516                         for (int i = 0; i < nThreads; i++) { model.collect(*mhps[i]); }
517                         model.finish();
518                 }
519
520                 // Relative error
521                 bChange = 0.0; totNum = 0;
522                 for (int i = 0; i <= M; i++)
523                         if (probv[i] >= 1e-7) {
524                                 change = fabs(theta[i] - probv[i]) / probv[i];
525                                 if (change >= STOP_CRITERIA) ++totNum;
526                                 if (bChange < change) bChange = change;
527                         }
528
529                 if (verbose) printf("ROUND = %d, SUM = %.15g, bChange = %f, totNum = %d\n", ROUND, sum, bChange, totNum);
530         } while (ROUND < MIN_ROUND || totNum > 0 && ROUND < MAX_ROUND);
531           //while (ROUND < MAX_ROUND);
532
533         if (totNum > 0) fprintf(stderr, "Warning: RSEM reaches %d iterations before meeting the convergence criteria.\n", MAX_ROUND);
534
535         //calculate expected effective lengths for each isoform
536         calcExpectedEffectiveLengths<ModelType>(model);
537
538         //correct theta vector
539         sum = theta[0];
540         for (int i = 1; i <= M; i++) 
541           if (eel[i] < EPSILON) { theta[i] = 0.0; }
542           else sum += theta[i];
543         if (sum < EPSILON) { fprintf(stderr, "No Expected Effective Length is no less than %.6g?!\n", MINEEL); exit(-1); }
544         for (int i = 0; i <= M; i++) theta[i] /= sum;
545
546         //generate output file used by Gibbs sampler
547         if (genGibbsOut) {
548                 if (model.getNeedCalcConPrb()) {
549                         for (int i = 0; i < nThreads; i++) {
550                                 rc = pthread_create(&threads[i], &attr, calcConProbs<ReadType, HitType, ModelType>, (void*)(&fparams[i]));
551                                 if (rc != 0) { fprintf(stderr, "Cannot create thread %d when generate files for Gibbs sampler! (numbered from 0)\n", i); exit(-1); }
552                         }
553                         for (int i = 0; i < nThreads; i++) {
554                                 rc = pthread_join(threads[i], &status);
555                                 if (rc != 0) { fprintf(stderr, "Cannot join thread %d when generate files for Gibbs sampler! (numbered from 0)\n", i); exit(-1); }
556                         }
557                 }
558                 model.setNeedCalcConPrb(false);
559
560                 sprintf(out_for_gibbs_F, "%s.ofg", imdName);
561                 FILE *fo = fopen(out_for_gibbs_F, "w");
562                 fprintf(fo, "%d %d\n", M, N0);
563                 for (int i = 0; i < nThreads; i++) {
564                         int numN = hitvs[i]->getN();
565                         for (int j = 0; j < numN; j++) {
566                                 int fr = hitvs[i]->getSAt(j);
567                                 int to = hitvs[i]->getSAt(j + 1);
568                                 int totNum = 0;
569
570                                 if (ncpvs[i][j] > 0.0) { ++totNum; fprintf(fo, "%d %.15g ", 0, ncpvs[i][j]); }
571                                 for (int k = fr; k < to; k++) {
572                                         HitType &hit = hitvs[i]->getHitAt(k);
573                                         if (hit.getConPrb() >= EPSILON) {
574                                                 ++totNum;
575                                                 fprintf(fo, "%d %.15g ", hit.getSid(), hit.getConPrb());
576                                         }
577                                 }
578
579                                 if (totNum > 0) { fprintf(fo, "\n"); }
580                         }
581                 }
582                 fclose(fo);
583         }
584
585         //calculate expected weights and counts using learned parameters
586         updateModel = false; calcExpectedWeights = true;
587         for (int i = 0; i < nThreads; i++) {
588                 rc = pthread_create(&threads[i], &attr, E_STEP<ReadType, HitType, ModelType>, (void*)(&fparams[i]));
589                 if (rc != 0) { fprintf(stderr, "Cannot create thread %d when calculate expected weights! (numbered from 0)\n", i); exit(-1); }
590                 //assert(rc == 0);
591         }
592         for (int i = 0; i < nThreads; i++) {
593                 rc = pthread_join(threads[i], &status);
594                 if (rc != 0) { fprintf(stderr, "Cannot join thread %d! (numbered from 0) when calculate expected weights!\n", i); exit(-1); }
595                 //assert(rc == 0);
596         }
597         model.setNeedCalcConPrb(false);
598         for (int i = 1; i < nThreads; i++) {
599                 for (int j = 0; j <= M; j++) {
600                         countvs[0][j] += countvs[i][j];
601                 }
602         }
603         countvs[0][0] += N0;
604
605         /* destroy attribute */
606         pthread_attr_destroy(&attr);
607
608         
609         //for all
610         double *mw = model.getMW();
611         sum = 0.0;
612         for (int i = 0; i <= M; i++) {
613           theta[i] = (mw[i] < EPSILON ? 0.0 : theta[i] / mw[i]);
614           sum += theta[i]; 
615         }
616         assert(sum >= EPSILON);
617         for (int i = 0; i <= M; i++) theta[i] /= sum;
618
619         writeResults<ModelType>(model, countvs[0]);
620
621         if (genBamF) {
622                 sprintf(outBamF, "%s.bam", outName);
623                 if (transcripts.getType() == 0) {
624                         sprintf(chr_list, "%s.chrlist", refName);
625                         pt_chr_list = (char*)(&chr_list);
626                 }
627
628                 BamWriter writer(inpSamType, inpSamF, pt_fn_list, outBamF, pt_chr_list);
629                 HitWrapper<HitType> wrapper(nThreads, hitvs);
630                 writer.work(wrapper, transcripts);
631         }
632
633         release<ReadType, HitType, ModelType>(readers, hitvs, ncpvs, mhps);
634 }
635
636 int main(int argc, char* argv[]) {
637         ifstream fin;
638         bool quiet = false;
639
640         if (argc < 5) {
641                 printf("Usage : rsem-run-em refName read_type imdName outName [-p #Threads] [-b samInpType samInpF has_fn_list_? [fn_list]] [-q] [--gibbs-out]\n\n");
642                 printf("  refName: reference name\n");
643                 printf("  read_type: 0 single read without quality score; 1 single read with quality score; 2 paired-end read without quality score; 3 paired-end read with quality score.\n");
644                 printf("  imdName: name for all upstream/downstream user-unseen files. (different files have different suffices)\n");
645                 printf("  outName: name for all output files. (different files have different suffices)\n");
646                 printf("  -p: number of threads which user wants to use. (default: 1)\n");
647                 printf("  -b: produce bam format output file. (default: off)\n");
648                 printf("  -q: set it quiet\n");
649                 printf("  --gibbs-out: generate output file use by Gibbs sampler. (default: off)\n");
650                 printf("// model parameters should be in imdName.mparams.\n");
651                 exit(-1);
652         }
653
654         time_t a = time(NULL);
655
656         strcpy(refName, argv[1]);
657         read_type = atoi(argv[2]);
658         strcpy(imdName, argv[3]);
659         strcpy(outName, argv[4]);
660
661         nThreads = 1;
662
663         genBamF = false;
664         genGibbsOut = false;
665         pt_fn_list = pt_chr_list = NULL;
666
667         for (int i = 5; i < argc; i++) {
668                 if (!strcmp(argv[i], "-p")) { nThreads = atoi(argv[i + 1]); }
669                 if (!strcmp(argv[i], "-b")) {
670                         genBamF = true;
671                         inpSamType = argv[i + 1][0];
672                         strcpy(inpSamF, argv[i + 2]);
673                         if (atoi(argv[i + 3]) == 1) {
674                                 strcpy(fn_list, argv[i + 4]);
675                                 pt_fn_list = (char*)(&fn_list);
676                         }
677                 }
678                 if (!strcmp(argv[i], "-q")) { quiet = true; }
679                 if (!strcmp(argv[i], "--gibbs-out")) { genGibbsOut = true; }
680         }
681         if (nThreads <= 0) { fprintf(stderr, "Number of threads should be bigger than 0!\n"); exit(-1); }
682         //assert(nThreads > 0);
683
684         verbose = !quiet;
685
686         //basic info loading
687         sprintf(refF, "%s.seq", refName);
688         refs.loadRefs(refF);
689         M = refs.getM();
690         sprintf(groupF, "%s.grp", refName);
691         gi.load(groupF);
692         m = gi.getm();
693
694         sprintf(tiF, "%s.ti", refName);
695         transcripts.readFrom(tiF);
696
697         sprintf(cntF, "%s.cnt", imdName);
698         fin.open(cntF);
699         if (!fin.is_open()) { fprintf(stderr, "Cannot open %s! It may not exist.\n", cntF); exit(-1); }
700         fin>>N0>>N1>>N2>>N_tot;
701         fin.close();
702
703         if (N1 <= 0) { fprintf(stderr, "There are no alignable reads!\n"); exit(-1); }
704
705         if (nThreads > N1) nThreads = N1;
706
707         //set model parameters
708         mparams.M = M;
709         mparams.N[0] = N0; mparams.N[1] = N1; mparams.N[2] = N2;
710         mparams.refs = &refs;
711
712         sprintf(mparamsF, "%s.mparams", imdName);
713         fin.open(mparamsF);
714         if (!fin.is_open()) { fprintf(stderr, "Cannot open %s! It may not exist.\n", mparamsF); exit(-1); }
715         fin>> mparams.minL>> mparams.maxL>> mparams.probF;
716         int val; // 0 or 1 , for estRSPD
717         fin>>val;
718         mparams.estRSPD = (val != 0);
719         fin>> mparams.B>> mparams.mate_minL>> mparams.mate_maxL>> mparams.mean>> mparams.sd;
720         fin>> mparams.seedLen;
721         fin.close();
722
723         //run EM
724         switch(read_type) {
725         case 0 : EM<SingleRead, SingleHit, SingleModel>(); break;
726         case 1 : EM<SingleReadQ, SingleHit, SingleQModel>(); break;
727         case 2 : EM<PairedEndRead, PairedEndHit, PairedEndModel>(); break;
728         case 3 : EM<PairedEndReadQ, PairedEndHit, PairedEndQModel>(); break;
729         default : fprintf(stderr, "Unknown Read Type!\n"); exit(-1);
730         }
731
732         time_t b = time(NULL);
733
734         printTimeUsed(a, b);
735
736         return 0;
737 }