]> git.donarmstrong.com Git - rsem.git/blob - EM.cpp
Normalized Read Fraction(nrf) are eliminated from outputs and posterior mean counts...
[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         //calculate tau values
319         double *tau = new double[M + 1];
320         memset(tau, 0, sizeof(double) * (M + 1));
321
322         denom = 0.0;
323         for (int i = 1; i <= M; i++) 
324           if (eel[i] >= EPSILON) {
325             tau[i] = theta[i] / eel[i];
326             denom += tau[i];
327           }   
328         if (denom <= 0) { fprintf(stderr, "No alignable reads?!\n"); exit(-1); }
329         //assert(denom > 0);
330         for (int i = 1; i <= M; i++) {
331                 tau[i] /= denom;
332         }
333
334         //isoform level results
335         sprintf(outF, "%s.iso_res", imdName);
336         fo = fopen(outF, "w");
337         for (int i = 1; i <= M; i++) {
338                 const Transcript& transcript = transcripts.getTranscriptAt(i);
339                 fprintf(fo, "%s%c", transcript.getTranscriptID().c_str(), (i < M ? '\t' : '\n'));
340         }
341         for (int i = 1; i <= M; i++)
342                 fprintf(fo, "%.2f%c", counts[i], (i < M ? '\t' : '\n'));
343         for (int i = 1; i <= M; i++)
344                 fprintf(fo, "%.15g%c", tau[i], (i < M ? '\t' : '\n'));
345         for (int i = 1; i <= M; i++) {
346                 const Transcript& transcript = transcripts.getTranscriptAt(i);
347                 fprintf(fo, "%s%c", transcript.getLeft().c_str(), (i < M ? '\t' : '\n'));
348         }
349         fclose(fo);
350
351         //gene level results
352         sprintf(outF, "%s.gene_res", imdName);
353         fo = fopen(outF, "w");
354         for (int i = 0; i < m; i++) {
355                 const string& gene_id = transcripts.getTranscriptAt(gi.spAt(i)).getGeneID();
356                 fprintf(fo, "%s%c", gene_id.c_str(), (i < m - 1 ? '\t' : '\n'));
357         }
358         for (int i = 0; i < m; i++) {
359                 double sumC = 0.0; // sum of counts
360                 int b = gi.spAt(i), e = gi.spAt(i + 1);
361                 for (int j = b; j < e; j++) sumC += counts[j];
362                 fprintf(fo, "%.2f%c", sumC, (i < m - 1 ? '\t' : '\n'));
363         }
364         for (int i = 0; i < m; i++) {
365                 double sumT = 0.0; // sum of tau values
366                 int b = gi.spAt(i), e = gi.spAt(i + 1);
367                 for (int j = b; j < e; j++) sumT += tau[j];
368                 fprintf(fo, "%.15g%c", sumT, (i < m - 1 ? '\t' : '\n'));
369         }
370         for (int i = 0; i < m; i++) {
371                 int b = gi.spAt(i), e = gi.spAt(i + 1);
372                 for (int j = b; j < e; j++) {
373                         fprintf(fo, "%s%c", transcripts.getTranscriptAt(j).getTranscriptID().c_str(), (j < e - 1 ? ',' : (i < m - 1 ? '\t' :'\n')));
374                 }
375         }
376         fclose(fo);
377
378         delete[] tau;
379
380         if (verbose) { printf("Expression Results are written!\n"); }
381 }
382
383 template<class ReadType, class HitType, class ModelType>
384 void release(ReadReader<ReadType> **readers, HitContainer<HitType> **hitvs, double **ncpvs, ModelType **mhps) {
385         delete[] probv;
386         for (int i = 0; i < nThreads; i++) {
387                 delete[] countvs[i];
388         }
389         delete[] countvs;
390
391         for (int i = 0; i < nThreads; i++) {
392                 delete readers[i];
393                 delete hitvs[i];
394                 delete[] ncpvs[i];
395                 delete mhps[i];
396         }
397         delete[] readers;
398         delete[] hitvs;
399         delete[] ncpvs;
400         delete[] mhps;
401 }
402
403 inline bool doesUpdateModel(int ROUND) {
404         //return false; // never update, for debugging only
405         return ROUND <= 20 || ROUND % 100 == 0;
406 }
407
408 //Including initialize, algorithm and results saving
409 template<class ReadType, class HitType, class ModelType>
410 void EM() {
411         int ROUND;
412         double sum;
413
414         double bChange = 0.0, change = 0.0; // bChange : biggest change
415         int totNum = 0;
416
417         ModelType model(mparams); //master model
418         ReadReader<ReadType> **readers;
419         HitContainer<HitType> **hitvs;
420         double **ncpvs;
421         ModelType **mhps; //model helpers
422
423         Params fparams[nThreads];
424         pthread_t threads[nThreads];
425         pthread_attr_t attr;
426         void *status;
427         int rc;
428
429         //initialize boolean variables
430         updateModel = calcExpectedWeights = false;
431
432         theta.clear();
433         theta.resize(M + 1, 0.0);
434         init<ReadType, HitType, ModelType>(readers, hitvs, ncpvs, mhps);
435
436         //set initial parameters
437         assert(N_tot > N2);
438         theta[0] = max(N0 * 1.0 / (N_tot - N2), 1e-8);
439         double val = (1.0 - theta[0]) / M;
440         for (int i = 1; i <= M; i++) theta[i] = val;
441
442         model.estimateFromReads(imdName);
443
444         for (int i = 0; i < nThreads; i++) {
445                 fparams[i].model = (void*)(&model);
446
447                 fparams[i].reader = (void*)readers[i];
448                 fparams[i].hitv = (void*)hitvs[i];
449                 fparams[i].ncpv = (void*)ncpvs[i];
450                 fparams[i].mhp = (void*)mhps[i];
451                 fparams[i].countv = (void*)countvs[i];
452         }
453
454         /* set thread attribute to be joinable */
455         pthread_attr_init(&attr);
456         pthread_attr_setdetachstate(&attr, PTHREAD_CREATE_JOINABLE);
457
458         ROUND = 0;
459         do {
460                 ++ROUND;
461
462                 updateModel = doesUpdateModel(ROUND);
463
464                 for (int i = 0; i <= M; i++) probv[i] = theta[i];
465
466                 //E step
467                 for (int i = 0; i < nThreads; i++) {
468                         rc = pthread_create(&threads[i], &attr, E_STEP<ReadType, HitType, ModelType>, (void*)(&fparams[i]));
469                         if (rc != 0) { fprintf(stderr, "Cannot create thread %d at ROUND %d! (numbered from 0)\n", i, ROUND); exit(-1); }
470                         //assert(rc == 0);
471                 }
472
473                 for (int i = 0; i < nThreads; i++) {
474                         rc = pthread_join(threads[i], &status);
475                         if (rc != 0) { fprintf(stderr, "Cannot join thread %d at ROUND %d! (numbered from 0)\n", i, ROUND); exit(-1); }
476                         //assert(rc == 0);
477                 }
478
479                 model.setNeedCalcConPrb(false);
480
481                 for (int i = 1; i < nThreads; i++) {
482                         for (int j = 0; j <= M; j++) {
483                                 countvs[0][j] += countvs[i][j];
484                         }
485                 }
486
487                 //add N0 noise reads
488                 countvs[0][0] += N0;
489
490                 //M step;
491                 sum = 0.0;
492                 for (int i = 0; i <= M; i++) sum += countvs[0][i];
493                 assert(sum >= EPSILON);
494                 for (int i = 0; i <= M; i++) theta[i] = countvs[0][i] / sum;
495
496                 if (updateModel) {
497                         model.init();
498                         for (int i = 0; i < nThreads; i++) { model.collect(*mhps[i]); }
499                         model.finish();
500                 }
501
502                 // Relative error
503                 bChange = 0.0; totNum = 0;
504                 for (int i = 0; i <= M; i++)
505                         if (probv[i] >= 1e-7) {
506                                 change = fabs(theta[i] - probv[i]) / probv[i];
507                                 if (change >= STOP_CRITERIA) ++totNum;
508                                 if (bChange < change) bChange = change;
509                         }
510
511                 if (verbose) printf("ROUND = %d, SUM = %.15g, bChange = %f, totNum = %d\n", ROUND, sum, bChange, totNum);
512         } while (ROUND < MIN_ROUND || totNum > 0 && ROUND < MAX_ROUND);
513           //while (ROUND < MAX_ROUND);
514
515         if (totNum > 0) fprintf(stderr, "Warning: RSEM reaches %d iterations before meeting the convergence criteria.\n", MAX_ROUND);
516
517         //calculate expected effective lengths for each isoform
518         calcExpectedEffectiveLengths<ModelType>(model);
519
520         //correct theta vector
521         sum = theta[0];
522         for (int i = 1; i <= M; i++) 
523           if (eel[i] < EPSILON) { theta[i] = 0.0; }
524           else sum += theta[i];
525         if (sum < EPSILON) { fprintf(stderr, "No Expected Effective Length is no less than %.6g?!\n", MINEEL); exit(-1); }
526         for (int i = 0; i <= M; i++) theta[i] /= sum;
527
528         //generate output file used by Gibbs sampler
529         if (genGibbsOut) {
530                 if (model.getNeedCalcConPrb()) {
531                         for (int i = 0; i < nThreads; i++) {
532                                 rc = pthread_create(&threads[i], &attr, calcConProbs<ReadType, HitType, ModelType>, (void*)(&fparams[i]));
533                                 if (rc != 0) { fprintf(stderr, "Cannot create thread %d when generate files for Gibbs sampler! (numbered from 0)\n", i); exit(-1); }
534                         }
535                         for (int i = 0; i < nThreads; i++) {
536                                 rc = pthread_join(threads[i], &status);
537                                 if (rc != 0) { fprintf(stderr, "Cannot join thread %d when generate files for Gibbs sampler! (numbered from 0)\n", i); exit(-1); }
538                         }
539                 }
540                 model.setNeedCalcConPrb(false);
541
542                 sprintf(out_for_gibbs_F, "%s.ofg", imdName);
543                 FILE *fo = fopen(out_for_gibbs_F, "w");
544                 fprintf(fo, "%d %d\n", M, N0);
545                 for (int i = 0; i < nThreads; i++) {
546                         int numN = hitvs[i]->getN();
547                         for (int j = 0; j < numN; j++) {
548                                 int fr = hitvs[i]->getSAt(j);
549                                 int to = hitvs[i]->getSAt(j + 1);
550                                 int totNum = 0;
551
552                                 if (ncpvs[i][j] > 0.0) { ++totNum; fprintf(fo, "%d %.15g ", 0, ncpvs[i][j]); }
553                                 for (int k = fr; k < to; k++) {
554                                         HitType &hit = hitvs[i]->getHitAt(k);
555                                         if (hit.getConPrb() >= EPSILON) {
556                                                 ++totNum;
557                                                 fprintf(fo, "%d %.15g ", hit.getSid(), hit.getConPrb());
558                                         }
559                                 }
560
561                                 if (totNum > 0) { fprintf(fo, "\n"); }
562                         }
563                 }
564                 fclose(fo);
565         }
566
567         //calculate expected weights and counts using learned parameters
568         updateModel = false; calcExpectedWeights = true;
569         for (int i = 0; i < nThreads; i++) {
570                 rc = pthread_create(&threads[i], &attr, E_STEP<ReadType, HitType, ModelType>, (void*)(&fparams[i]));
571                 if (rc != 0) { fprintf(stderr, "Cannot create thread %d when calculate expected weights! (numbered from 0)\n", i); exit(-1); }
572                 //assert(rc == 0);
573         }
574         for (int i = 0; i < nThreads; i++) {
575                 rc = pthread_join(threads[i], &status);
576                 if (rc != 0) { fprintf(stderr, "Cannot join thread %d! (numbered from 0) when calculate expected weights!\n", i); exit(-1); }
577                 //assert(rc == 0);
578         }
579         model.setNeedCalcConPrb(false);
580         for (int i = 1; i < nThreads; i++) {
581                 for (int j = 0; j <= M; j++) {
582                         countvs[0][j] += countvs[i][j];
583                 }
584         }
585         countvs[0][0] += N0;
586
587         /* destroy attribute */
588         pthread_attr_destroy(&attr);
589
590         
591         //for all
592         double *mw = model.getMW();
593         sum = 0.0;
594         for (int i = 0; i <= M; i++) {
595           theta[i] = (mw[i] < EPSILON ? 0.0 : theta[i] / mw[i]);
596           sum += theta[i]; 
597         }
598         assert(sum >= EPSILON);
599         for (int i = 0; i <= M; i++) theta[i] /= sum;
600
601         writeResults<ModelType>(model, countvs[0]);
602
603         if (genBamF) {
604                 sprintf(outBamF, "%s.bam", outName);
605                 if (transcripts.getType() == 0) {
606                         sprintf(chr_list, "%s.chrlist", refName);
607                         pt_chr_list = (char*)(&chr_list);
608                 }
609
610                 BamWriter writer(inpSamType, inpSamF, pt_fn_list, outBamF, pt_chr_list);
611                 HitWrapper<HitType> wrapper(nThreads, hitvs);
612                 writer.work(wrapper, transcripts);
613         }
614
615         release<ReadType, HitType, ModelType>(readers, hitvs, ncpvs, mhps);
616 }
617
618 int main(int argc, char* argv[]) {
619         ifstream fin;
620         bool quiet = false;
621
622         if (argc < 5) {
623                 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");
624                 printf("  refName: reference name\n");
625                 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");
626                 printf("  imdName: name for all upstream/downstream user-unseen files. (different files have different suffices)\n");
627                 printf("  outName: name for all output files. (different files have different suffices)\n");
628                 printf("  -p: number of threads which user wants to use. (default: 1)\n");
629                 printf("  -b: produce bam format output file. (default: off)\n");
630                 printf("  -q: set it quiet\n");
631                 printf("  --gibbs-out: generate output file use by Gibbs sampler. (default: off)\n");
632                 printf("// model parameters should be in imdName.mparams.\n");
633                 exit(-1);
634         }
635
636         time_t a = time(NULL);
637
638         strcpy(refName, argv[1]);
639         read_type = atoi(argv[2]);
640         strcpy(imdName, argv[3]);
641         strcpy(outName, argv[4]);
642
643         nThreads = 1;
644
645         genBamF = false;
646         genGibbsOut = false;
647         pt_fn_list = pt_chr_list = NULL;
648
649         for (int i = 5; i < argc; i++) {
650                 if (!strcmp(argv[i], "-p")) { nThreads = atoi(argv[i + 1]); }
651                 if (!strcmp(argv[i], "-b")) {
652                         genBamF = true;
653                         inpSamType = argv[i + 1][0];
654                         strcpy(inpSamF, argv[i + 2]);
655                         if (atoi(argv[i + 3]) == 1) {
656                                 strcpy(fn_list, argv[i + 4]);
657                                 pt_fn_list = (char*)(&fn_list);
658                         }
659                 }
660                 if (!strcmp(argv[i], "-q")) { quiet = true; }
661                 if (!strcmp(argv[i], "--gibbs-out")) { genGibbsOut = true; }
662         }
663         if (nThreads <= 0) { fprintf(stderr, "Number of threads should be bigger than 0!\n"); exit(-1); }
664         //assert(nThreads > 0);
665
666         verbose = !quiet;
667
668         //basic info loading
669         sprintf(refF, "%s.seq", refName);
670         refs.loadRefs(refF);
671         M = refs.getM();
672         sprintf(groupF, "%s.grp", refName);
673         gi.load(groupF);
674         m = gi.getm();
675
676         sprintf(tiF, "%s.ti", refName);
677         transcripts.readFrom(tiF);
678
679         sprintf(cntF, "%s.cnt", imdName);
680         fin.open(cntF);
681         if (!fin.is_open()) { fprintf(stderr, "Cannot open %s! It may not exist.\n", cntF); exit(-1); }
682         fin>>N0>>N1>>N2>>N_tot;
683         fin.close();
684
685         if (N1 <= 0) { fprintf(stderr, "There are no alignable reads!\n"); exit(-1); }
686
687         if (nThreads > N1) nThreads = N1;
688
689         //set model parameters
690         mparams.M = M;
691         mparams.N[0] = N0; mparams.N[1] = N1; mparams.N[2] = N2;
692         mparams.refs = &refs;
693
694         sprintf(mparamsF, "%s.mparams", imdName);
695         fin.open(mparamsF);
696         if (!fin.is_open()) { fprintf(stderr, "Cannot open %s! It may not exist.\n", mparamsF); exit(-1); }
697         fin>> mparams.minL>> mparams.maxL>> mparams.probF;
698         int val; // 0 or 1 , for estRSPD
699         fin>>val;
700         mparams.estRSPD = (val != 0);
701         fin>> mparams.B>> mparams.mate_minL>> mparams.mate_maxL>> mparams.mean>> mparams.sd;
702         fin>> mparams.seedLen;
703         fin.close();
704
705         //run EM
706         switch(read_type) {
707         case 0 : EM<SingleRead, SingleHit, SingleModel>(); break;
708         case 1 : EM<SingleReadQ, SingleHit, SingleQModel>(); break;
709         case 2 : EM<PairedEndRead, PairedEndHit, PairedEndModel>(); break;
710         case 3 : EM<PairedEndReadQ, PairedEndHit, PairedEndQModel>(); break;
711         default : fprintf(stderr, "Unknown Read Type!\n"); exit(-1);
712         }
713
714         time_t b = time(NULL);
715
716         printTimeUsed(a, b);
717
718         return 0;
719 }