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