-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathtask_4_5.C
499 lines (402 loc) · 18.5 KB
/
task_4_5.C
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
#include "common.h"
// build peak backgrounds
void build_pdf_peak(RooWorkspace *wspace, uint cate, double bdt_min)
{
gROOT->SetBatch();
using namespace RooFit;
vector<TString> decay;
vector<double> yield, yield_err;
// B0 -> mu mu
decay.push_back("bdmmMc");
yield.push_back(effyield[cate][N_bdmmMc]);
yield_err.push_back(effyield[cate][dN_bdmmMc]);
// Bs -> hadron hadron (Bs-> kk ; Bs-> k pi; Bs-> pi pi)
decay.push_back("bstohhMcBg");
yield.push_back(effyield[cate][N_bstohhMcBg]);
yield_err.push_back(effyield[cate][dN_bstohhMcBg]);
// B -> hadron hadron (B-> kk ; B-> k pi; B-> pi pi)
decay.push_back("bdtohhMcBg");
yield.push_back(effyield[cate][N_bdtohhMcBg]);
yield_err.push_back(effyield[cate][dN_bdtohhMcBg]);
RooRealVar m("m","",4.9,5.9);
RooRealVar wgt("wgt","",1.,0.,1000.);
RooDataSet *rds = new RooDataSet("rds","",RooArgSet(m,wgt),"wgt");
double sum_weight = 0.;
double sum_weight_err = 0.;
rds->weightError(RooAbsData::SumW2);
double weight_max = 0.; // get the largest weight value across samples
for(int proc=0; proc<(int)decay.size(); proc++) {
TFile *fin = new TFile("/eos/uscms/store/user/cmsdas/2025/long_exercises/long-ex-bs-mumu/" + decay[proc]+".root");
TTree *tin = (TTree*)fin->Get(decay[proc]);
double weight = yield[proc]/(double)tin->GetEntries(Form("cate==%d",cate));
if (weight>weight_max) weight_max = weight;
delete fin;
}
for(int proc=0; proc<(int)decay.size(); proc++) {
TFile *fin = new TFile("/eos/uscms/store/user/cmsdas/2025/long_exercises/long-ex-bs-mumu/" + decay[proc]+".root");
TTree *tin = (TTree*)fin->Get(decay[proc]);
unsigned int cate_t;
float m_t,bdt_t;
tin->SetBranchAddress("cate",&cate_t);
tin->SetBranchAddress("m",&m_t);
tin->SetBranchAddress("bdt",&bdt_t);
double weight = yield[proc]/(double)tin->GetEntries(Form("cate==%d",cate));
double weight_err = yield_err[proc]/(double)tin->GetEntries(Form("cate==%d",cate))/weight_max;
std::cout<<"weight "<<weight<<" weight_err"<<weight_err<<std::endl;
for(int evt=0; evt<tin->GetEntries(); evt++) {
tin->GetEntry(evt);
if (cate_t!=cate) continue;
if (bdt_t<=bdt_min) continue;
m.setVal(m_t);
wgt.setVal(weight/weight_max); // rescale the event with largest weight to be 1
rds->add(RooArgSet(m,wgt),weight/weight_max);
sum_weight += weight;
sum_weight_err += weight_err; // systematics; linear sum
}
delete fin;
}
RooKeysPdf pdf(Form("pdf_peak_%d",cate), "", m, *rds, RooKeysPdf::NoMirror, 2.0);
RooPlot *frame = m.frame(Title(" "),Bins(100));
rds->plotOn(frame, Name("t_rds"));
pdf.plotOn(frame, Name("t_pdf"), LineWidth(3));
TCanvas* canvas = new TCanvas("canvas", "", 600, 600);
canvas->SetMargin(0.15,0.06,0.13,0.07);
frame->GetYaxis()->SetTitleOffset(1.50);
frame->GetYaxis()->SetTitle("Entries / 0.01 GeV");
frame->GetXaxis()->SetTitleOffset(1.15);
frame->GetXaxis()->SetLabelOffset(0.01);
frame->GetXaxis()->SetTitle("M(#mu#mu) [GeV]");
frame->GetXaxis()->SetTitleSize(0.043);
frame->GetYaxis()->SetTitleSize(0.043);
frame->Draw();
TLegend* leg = new TLegend(0.58,0.77,0.93,0.92);
leg->SetFillStyle(0);
leg->SetLineWidth(0);
leg->SetHeader(Form("Category %d",cate));
leg->AddEntry(frame->findObject("t_rds"),"Simluation","EP");
leg->AddEntry(frame->findObject("t_pdf"),"PDF","L");
leg->Draw();
canvas->Print(Form("task_4_5a_%d.pdf",cate));
cout << "Category: " << cate << endl;
cout << "BDT min: " << bdt_min << endl;
cout << "Sum of weights: " << sum_weight << " +- " << sum_weight_err << endl;
wspace->import(pdf);
wspace->import(RooRealVar(Form("N_peak_%d",cate),"",sum_weight));
delete frame;
delete rds;
delete leg;
delete canvas;
}
// build semileptonic backgrounds
void build_pdf_semi(RooWorkspace *wspace, uint cate, double bdt_min)
{
vector<TString> decay;
vector<double> yield, yield_err;
decay.push_back("bskmunuMcBg");
yield.push_back(effyield[cate][N_bskmunuMcBg]);
yield_err.push_back(effyield[cate][dN_bskmunuMcBg]);
decay.push_back("bdpimunuMcBg");
yield.push_back(effyield[cate][N_bdpimunuMcBg]);
yield_err.push_back(effyield[cate][dN_bdpimunuMcBg]);
decay.push_back("bdpimumuMcBg");
yield.push_back(effyield[cate][N_bdpimumuMcBg]);
yield_err.push_back(effyield[cate][dN_bdpimumuMcBg]);
decay.push_back("bupimumuMcBg");
yield.push_back(effyield[cate][N_bupimumuMcBg]);
yield_err.push_back(effyield[cate][dN_bupimumuMcBg]);
RooRealVar m("m","",4.9,5.9);
RooRealVar wgt("wgt","",1.,0.,1000.);
RooDataSet *rds = new RooDataSet("rds","",RooArgSet(m,wgt),"wgt");
double sum_weight = 0.;
double sum_weight_err = 0.;
for(int proc=0; proc<(int)decay.size(); proc++) {
TFile *fin = new TFile("/eos/uscms/store/user/cmsdas/2025/long_exercises/long-ex-bs-mumu/"+decay[proc]+".root");
TTree *tin = (TTree*)fin->Get(decay[proc]);
uint cate_t;
float m_t,bdt_t;
tin->SetBranchAddress("cate",&cate_t);
tin->SetBranchAddress("m",&m_t);
tin->SetBranchAddress("bdt",&bdt_t);
double weight = yield[proc]/(double)tin->GetEntries(Form("cate==%d",cate));
double weight_err = yield_err[proc]/(double)tin->GetEntries(Form("cate==%d",cate));
for(int evt=0; evt<tin->GetEntries(); evt++) {
tin->GetEntry(evt);
if (cate_t!=cate) continue;
if (bdt_t<=bdt_min) continue;
m.setVal(m_t);
wgt.setVal(weight);
rds->add(RooArgSet(m,wgt),weight);
sum_weight += weight;
sum_weight_err += weight_err; // systematics; linear sum
}
delete fin;
}
RooKeysPdf pdf(Form("pdf_semi_%d",cate), "", m, *rds, RooKeysPdf::MirrorLeft, 2.0);
RooPlot *frame = m.frame(Title(" "),Bins(50));
rds->plotOn(frame, Name("t_rds"));
pdf.plotOn(frame, Name("t_pdf"), LineWidth(3));
TCanvas* canvas = new TCanvas("canvas", "", 600, 600);
canvas->SetMargin(0.15,0.06,0.13,0.07);
frame->GetYaxis()->SetTitleOffset(1.50);
frame->GetYaxis()->SetTitle("Entries / 0.02 GeV");
frame->GetXaxis()->SetTitleOffset(1.15);
frame->GetXaxis()->SetLabelOffset(0.01);
frame->GetXaxis()->SetTitle("M(#mu#mu) [GeV]");
frame->GetXaxis()->SetTitleSize(0.043);
frame->GetYaxis()->SetTitleSize(0.043);
frame->Draw();
TLegend* leg = new TLegend(0.58,0.77,0.93,0.92);
leg->SetFillStyle(0);
leg->SetLineWidth(0);
leg->SetHeader(Form("Category %d",cate));
leg->AddEntry(frame->findObject("t_rds"),"Simluation","EP");
leg->AddEntry(frame->findObject("t_pdf"),"PDF","L");
leg->Draw();
canvas->Print(Form("task_4_5b_%d.pdf",cate));
cout << "Category: " << cate << endl;
cout << "BDT min: " << bdt_min << endl;
cout << "Sum of weights: " << sum_weight << " +- " << sum_weight_err << endl;
wspace->import(pdf);
wspace->import(RooRealVar(Form("N_semi_%d",cate),"",sum_weight));
delete frame;
delete rds;
delete leg;
delete canvas;
}
// build signal PDF, fit it to data, and fix the parameters
void build_pdf_bs(RooWorkspace *wspace, uint cate, double bdt_min)
{
RooRealVar m("m","",4.9,5.9);
RooDataSet *rds = new RooDataSet("rds","",RooArgSet(m));
TFile *fin = new TFile("/eos/uscms/store/user/cmsdas/2025/long_exercises/long-ex-bs-mumu/bsmmMc.root");
TTree *tin = (TTree*)fin->Get("bsmmMc");
uint cate_t;
float m_t,bdt_t;
tin->SetBranchAddress("cate",&cate_t);
tin->SetBranchAddress("m",&m_t);
tin->SetBranchAddress("bdt",&bdt_t);
int evtcount[2] = {0,0};
for(int evt=0; evt<tin->GetEntries(); evt++) {
tin->GetEntry(evt);
if (cate_t!=cate) continue;
evtcount[0]++;
if (bdt_t<=bdt_min) continue;
evtcount[1]++;
m.setVal(m_t);
rds->add(RooArgSet(m));
}
double eff = (double)evtcount[1]/(double)evtcount[0] * effyield[cate][Eff_bsmmMc];
double eff_error = effyield[cate][dEff_bsmmMc]/effyield[cate][Eff_bsmmMc]*eff; // ignore MC statistical error
delete fin;
RooRealVar bs_mean1(Form("bs_%d_mean1",cate),"",5.37,5.2,5.5);
RooRealVar bs_mean2(Form("bs_%d_mean2",cate),"",5.37,5.2,5.5);
RooRealVar bs_sigma1(Form("bs_%d_sigma1",cate),"",0.030,0.005,0.060);
RooRealVar bs_sigma2(Form("bs_%d_sigma2",cate),"",0.080,0.040,0.200);
RooRealVar bs_cbalpha(Form("bs_%d_cbalpha",cate),"",1.,0.,4.);
RooRealVar bs_cbn(Form("bs_%d_cbn",cate),"",1.,0.,4.);
RooRealVar bs_frac(Form("bs_%d_frac",cate),"",0.7,0.,1.);
RooGaussian bs_gaus(Form("bs_%d_gaus",cate), "", m, bs_mean1, bs_sigma1);
RooCBShape bs_cbline(Form("bs_%d_cbline",cate), "", m, bs_mean2, bs_sigma2, bs_cbalpha, bs_cbn);
RooAddPdf pdf(Form("pdf_bs_%d",cate),"",RooArgList(bs_gaus,bs_cbline),RooArgList(bs_frac));
pdf.fitTo(*rds);
RooPlot *frame = m.frame(Title(" "),Bins(100));
rds->plotOn(frame, Name("t_rds"));
pdf.plotOn(frame, Name("t_pdf"), LineWidth(3));
TCanvas* canvas = new TCanvas("canvas", "", 600, 600);
canvas->SetMargin(0.15,0.06,0.13,0.07);
frame->GetYaxis()->SetTitleOffset(1.50);
frame->GetYaxis()->SetTitle("Entries / 0.01 GeV");
frame->GetXaxis()->SetTitleOffset(1.15);
frame->GetXaxis()->SetLabelOffset(0.01);
frame->GetXaxis()->SetTitle("M(#mu#mu) [GeV]");
frame->GetXaxis()->SetTitleSize(0.043);
frame->GetYaxis()->SetTitleSize(0.043);
frame->Draw();
TLegend* leg = new TLegend(0.58,0.77,0.93,0.92);
leg->SetFillStyle(0);
leg->SetLineWidth(0);
leg->SetHeader(Form("Category %d",cate));
leg->AddEntry(frame->findObject("t_rds"),"Simluation","EP");
leg->AddEntry(frame->findObject("t_pdf"),"PDF","L");
leg->Draw();
canvas->Print(Form("task_4_5c_%d.pdf",cate));
cout << "Category: " << cate << endl;
cout << "BDT min: " << bdt_min << endl;
bs_mean1.setConstant(true);
bs_mean2.setConstant(true);
bs_sigma1.setConstant(true);
bs_sigma2.setConstant(true);
bs_cbalpha.setConstant(true);
bs_cbn.setConstant(true);
bs_frac.setConstant(true);
wspace->import(pdf);
wspace->import(RooRealVar(Form("Eff_bs_%d",cate),"",eff));
delete frame;
delete rds;
delete leg;
delete canvas;
}
//normalisation channel
void fit_bupsik(RooWorkspace *wspace, uint cate, double bdt_min)
{
RooRealVar m("m","",5.0,5.8);
RooRealVar wgt("wgt","",1.,0.,1000.);
RooDataSet *rds_data = new RooDataSet("rds_data","",RooArgSet(m,wgt),"wgt");
RooDataSet *rds_mc = new RooDataSet("rds_mc","",RooArgSet(m));
TFile *fin = new TFile("/eos/uscms/store/user/cmsdas/2025/long_exercises/long-ex-bs-mumu/bupsikData.root");
TTree *tin = (TTree*)fin->Get("bupsikData");
uint cate_t;
float wgt_t,m_t,bdt_t;
tin->SetBranchAddress("cate",&cate_t);
tin->SetBranchAddress("wgt",&wgt_t);
tin->SetBranchAddress("m",&m_t);
for(int evt=0; evt<tin->GetEntries(); evt++) {
tin->GetEntry(evt);
if (cate_t!=cate) continue;
if (m_t<5.0 || m_t>=5.8) continue;
m.setVal(m_t);
wgt.setVal(wgt_t);
rds_data->add(RooArgSet(m,wgt),wgt_t);
}
delete fin;
fin = new TFile("/eos/uscms/store/user/cmsdas/2025/long_exercises/long-ex-bs-mumu/bupsikMc.root");
tin = (TTree*)fin->Get("bupsikMc");
tin->SetBranchAddress("cate",&cate_t);
tin->SetBranchAddress("m",&m_t);
int evtcount[2] = {0,0};
for(int evt=0; evt<tin->GetEntries(); evt++) {
tin->GetEntry(evt);
if (cate_t!=cate) continue;
evtcount[0]++;
if (m_t<5.0 || m_t>=5.8) continue;
evtcount[1]++;
m.setVal(m_t);
rds_mc->add(RooArgSet(m));
}
double eff = (double)evtcount[1]/(double)evtcount[0] * effyield[cate][Eff_bupsikMc];
double eff_error = effyield[cate][dEff_bupsikMc]/effyield[cate][Eff_bupsikMc]*eff; // ignore MC statistics error
delete fin;
RooRealVar sigmc_mean1("sigmc_mean1","",5.28,5.2,5.4);
RooRealVar sigmc_mean2("sigmc_mean2","",5.28,5.2,5.4);
RooRealVar sigmc_sigma1("sigmc_sigma1","",0.030,0.005,0.060);
RooRealVar sigmc_sigma2("sigmc_sigma2","",0.080,0.040,0.200);
RooRealVar sig_frac("sig_frac","",0.9,0.5,1.0);
RooGaussian sigmc_g1("sig_g1","",m,sigmc_mean1,sigmc_sigma1);
RooGaussian sigmc_g2("sig_g2","",m,sigmc_mean2,sigmc_sigma2);
RooAddPdf pdf_sigmc("pdf_sigmc","",RooArgList(sigmc_g1,sigmc_g2),RooArgList(sig_frac));
pdf_sigmc.fitTo(*rds_mc);
RooPlot *frame1 = m.frame(Title(" "),Bins(80));
rds_mc->plotOn(frame1, Name("t_rds_mc"));
pdf_sigmc.plotOn(frame1, Name("t_pdf_sigmc"), LineWidth(3));
TCanvas* canvas1 = new TCanvas("canvas1", "", 600, 600);
canvas1->SetMargin(0.15,0.06,0.13,0.07);
frame1->GetYaxis()->SetTitleOffset(1.50);
frame1->GetYaxis()->SetTitle("Entries / 0.01 GeV");
frame1->GetXaxis()->SetTitleOffset(1.15);
frame1->GetXaxis()->SetLabelOffset(0.01);
frame1->GetXaxis()->SetTitle("M(#mu#muK) [GeV]");
frame1->GetXaxis()->SetTitleSize(0.043);
frame1->GetYaxis()->SetTitleSize(0.043);
frame1->Draw();
TLegend* leg1 = new TLegend(0.58,0.77,0.93,0.92);
leg1->SetFillStyle(0);
leg1->SetLineWidth(0);
leg1->SetHeader(Form("Category %d",cate));
leg1->AddEntry(frame1->findObject("t_rds_mc"),"Simluation","EP");
leg1->AddEntry(frame1->findObject("t_pdf_sigmc"),"Fit","L");
leg1->Draw();
canvas1->Print(Form("task_4_5d_%d.pdf",cate));
sigmc_mean1.setConstant(true);
sigmc_mean2.setConstant(true);
sigmc_sigma1.setConstant(true);
sigmc_sigma2.setConstant(true);
sig_frac.setConstant(true);
RooRealVar sig_shift("sig_shift","",0.,-0.02,0.02);
RooRealVar sig_scale("sig_scale","",1.,0.8,1.2);
RooAddition sig_mean1("sig_mean1","",RooArgList(sigmc_mean1,sig_shift));
RooAddition sig_mean2("sig_mean2","",RooArgList(sigmc_mean2,sig_shift));
RooProduct sig_sigma1("sig_sigma1","",RooArgList(sigmc_sigma1,sig_scale));
RooProduct sig_sigma2("sig_sigma2","",RooArgList(sigmc_sigma2,sig_scale));
RooGaussian sig_g1("sig_g1","",m,sig_mean1,sig_sigma1);
RooGaussian sig_g2("sig_g2","",m,sig_mean2,sig_sigma2);
RooAddPdf pdf_sig("pdf_sig","",RooArgList(sig_g1,sig_g2),RooArgList(sig_frac));
RooRealVar comb_coeff("comb_coeff","",-1.2,-10.,10.);
RooExponential pdf_comb("pdf_comb","",m,comb_coeff);
RooRealVar jpsix_scale("jpsix_scale","",0.02,0.001,0.08);
RooRealVar jpsix_shift("jpsix_shift","",5.13,5.12,5.16);
RooGenericPdf pdf_jpsix("pdf_jpsix","","TMath::Erfc((@0-@1)/@2)",RooArgList(m,jpsix_shift,jpsix_scale));
double n_comb_guess = rds_data->sumEntries("m>5.4")*2.;
double n_sig_guess = rds_data->sumEntries("m>5.18&&m<5.38")-n_comb_guess/4.;
double n_jpsix_guess = rds_data->sumEntries("m<5.18")-n_comb_guess*0.18/0.8;
RooRealVar n_sig("n_sig","",n_sig_guess,0.,rds_data->sumEntries());
RooRealVar n_comb("n_comb","",n_comb_guess,0.,rds_data->sumEntries());
RooRealVar n_jpsix("n_jpsix","",n_jpsix_guess,0.,rds_data->sumEntries());
RooAddPdf model("model","",RooArgList(pdf_sig,pdf_comb,pdf_jpsix),RooArgList(n_sig,n_comb,n_jpsix));
model.fitTo(*rds_data, SumW2Error(true));
RooPlot *frame2 = m.frame(Title(" "),Bins(80));
rds_data->plotOn(frame2, Name("t_rds_data"));
model.plotOn(frame2, Name("t_model"), LineWidth(3));
model.plotOn(frame2, Name("t_pdf_comb"), Components("pdf_comb"), LineWidth(3), LineStyle(2), LineColor(kGray+1));
TCanvas* canvas2 = new TCanvas("canvas2", "", 600, 600);
canvas2->SetMargin(0.15,0.06,0.13,0.07);
frame2->GetYaxis()->SetTitleOffset(1.50);
frame2->GetYaxis()->SetTitle("Entries / 0.01 GeV");
frame2->GetXaxis()->SetTitleOffset(1.15);
frame2->GetXaxis()->SetLabelOffset(0.01);
frame2->GetXaxis()->SetTitle("M(#mu#muK) [GeV]");
frame2->GetXaxis()->SetTitleSize(0.043);
frame2->GetYaxis()->SetTitleSize(0.043);
frame2->Draw();
TLegend* leg2 = new TLegend(0.58,0.77,0.93,0.92);
leg2->SetFillStyle(0);
leg2->SetLineWidth(0);
leg2->SetHeader(Form("Category %d",cate));
leg2->AddEntry(frame2->findObject("t_rds_data"),"Data","EP");
leg2->AddEntry(frame2->findObject("t_model"),"Fit","L");
leg2->AddEntry(frame2->findObject("t_pdf_comb"),"Combinatorial bkg.","L");
leg2->Draw();
canvas2->Print(Form("task_4_5e_%d.pdf",cate));
cout << "Category: " << cate << endl;
cout << "Selection efficiency: " << eff << " +- " << eff_error << endl;
cout << "Observed yield: " << n_sig.getVal() << " +- " << n_sig.getError() << endl;
wspace->import(RooRealVar(Form("N_bu_%d",cate),"",n_sig.getVal()));
wspace->import(RooRealVar(Form("Eff_bu_%d",cate),"",eff));
delete frame1;
delete frame2;
delete rds_mc;
delete rds_data;
delete leg1;
delete leg2;
delete canvas1;
delete canvas2;
}
// combinatorial bkg PDF model
void build_pdf_comb(RooWorkspace *wspace, uint cate, double bdt_min)
{
RooRealVar m("m","",4.9,5.9);
RooRealVar comb_B1(Form("comb_B1_%d",cate), "", 0.5, 0. , 1);
RooFormulaVar comb_B2(Form("comb_B2_%d",cate), "", "1.-@0", RooArgList(comb_B1));
RooBernstein pdf(Form("pdf_comb_%d",cate), "", m, RooArgList(comb_B1, comb_B2));
TFile *fin = new TFile("/eos/uscms/store/user/cmsdas/2025/long_exercises/long-ex-bs-mumu/bmmData-blind.root");
TTree *tin = (TTree*)fin->Get("bmmData");
double N_comb_guess = (double)tin->GetEntries(Form("cate==%d&&bdt>%g&&m>5.45",cate,bdt_min));
N_comb_guess *= 1.0/0.45; // scale to full mass region
delete fin;
cout << "Category: " << cate << endl;
wspace->import(pdf);
wspace->import(RooRealVar(Form("N_comb_%d",cate),"",N_comb_guess,0.,N_comb_guess*10.));
}
void task_4_5()
{
RooWorkspace *wspace = new RooWorkspace("wspace");
for(uint cate=0; cate<N_Categories; cate++) {
double bdt_min = 0.80;
build_pdf_peak(wspace,cate,bdt_min);
build_pdf_semi(wspace,cate,bdt_min);
build_pdf_bs(wspace,cate,bdt_min);
build_pdf_comb(wspace,cate,bdt_min);
fit_bupsik(wspace, cate, bdt_min);
}
wspace->writeToFile("wspace.root");
delete wspace;
}