00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include <TMath.h>
00011 #include <TCanvas.h>
00012 #include <TRandom3.h>
00013 #include <TFitter.h>
00014 #include <TF1.h>
00015 #include <TStyle.h>
00016 #include <TVector.h>
00017 #include <TGraph.h>
00018
00019 #include <TUnfoldSys.h>
00020
00021 using namespace std;
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045 TRandom *rnd=0;
00046
00047 Double_t GenerateEvent(const Double_t *parm,
00048 const Double_t *triggerParm,
00049 Double_t *intLumi,
00050 Bool_t *triggerFlag,
00051 Double_t *ptGen,Int_t *iType)
00052 {
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093 Double_t ptObs;
00094 Bool_t isTriggered=kFALSE;
00095 do {
00096 Int_t itype;
00097 Double_t ptgen;
00098 Double_t f=rnd->Rndm();
00099
00100 itype=0;
00101 if(f<parm[0]) itype=1;
00102 else if(f<parm[0]+parm[1]) itype=2;
00103
00104
00105 Double_t a=parm[4+itype];
00106 Double_t a1=a+1.0;
00107 Double_t t=rnd->Rndm();
00108 if(a1 == 0.0) {
00109 Double_t x0=TMath::Log(parm[2]);
00110 ptgen=TMath::Exp(t*(TMath::Log(parm[3])-x0)+x0);
00111 } else {
00112 Double_t x0=pow(parm[2],a1);
00113 ptgen=pow(t*(pow(parm[3],a1)-x0)+x0,1./a1);
00114 }
00115 if(iType) *iType=itype;
00116 if(ptGen) *ptGen=ptgen;
00117
00118
00119 Double_t sigma=
00120 TMath::Sqrt(parm[7]*parm[7]*ptgen+parm[8]*parm[8]*ptgen*ptgen);
00121 ptObs=rnd->BreitWigner(ptgen,sigma);
00122
00123
00124 Double_t triggerProb =
00125 triggerParm[2]/(1.+TMath::Exp((triggerParm[0]-ptObs)/triggerParm[1]));
00126 isTriggered= rnd->Rndm()<triggerProb;
00127 (*intLumi) ++;
00128 } while((!triggerFlag) && (!isTriggered));
00129
00130 if(triggerFlag) *triggerFlag=isTriggered;
00131 return ptObs;
00132 }
00133
00134
00135 void testUnfold3()
00136 {
00137
00138 TH1::SetDefaultSumw2();
00139
00140
00141 gStyle->SetOptFit(1111);
00142
00143
00144 rnd=new TRandom3();
00145
00146
00147 Double_t const lumiData= 30000;
00148 Double_t const lumiMC =1000000;
00149
00150
00151
00152 Int_t const nDet=24;
00153 Double_t const xminDet=4.0;
00154 Double_t const xmaxDet=28.0;
00155
00156 Int_t const nGen=10;
00157 Double_t const xminGen= 6.0;
00158 Double_t const xmaxGen=26.0;
00159
00160
00161
00162
00163
00164
00165
00166
00167 static Double_t parm_data[]={
00168 0.05,
00169 0.05,
00170 3.5,
00171 100.,
00172 -3.0,
00173 0.1,
00174 -0.5,
00175 0.2,
00176 0.01,
00177 };
00178 static Double_t triggerParm_data[]={8.0,
00179 4.0,
00180 0.95
00181 };
00182
00183
00184
00185 TH1D *histDetDATA=new TH1D("PT(data)",";Pt(obs)",nDet,xminDet,xmaxDet);
00186
00187
00188 TH1D *histDetDATAbbb=new TH1D("PT(data,bbb)",";Pt(obs,bbb)",
00189 nGen,xminGen,xmaxGen);
00190
00191 TH1D *histGenDATA=new TH1D("PT(MC,signal,gen)",";Pt(gen)",
00192 nGen,xminGen,xmaxGen);
00193
00194 Double_t intLumi=0.0;
00195 while(intLumi<lumiData) {
00196 Int_t iTypeGen;
00197 Bool_t isTriggered;
00198 Double_t ptGen;
00199 Double_t ptObs=GenerateEvent(parm_data,triggerParm_data,&intLumi,
00200 &isTriggered,&ptGen,&iTypeGen);
00201 if(isTriggered) {
00202 histDetDATA->Fill(ptObs);
00203 histDetDATAbbb->Fill(ptObs);
00204 }
00205
00206
00207 if(iTypeGen==0) histGenDATA->Fill(ptGen);
00208 }
00209
00210
00211
00212
00213 static Double_t parm_MC[]={
00214 0.05,
00215 0.05,
00216 3.5,
00217 100.,
00218 -4.0,
00219
00220 0.1,
00221 -0.5,
00222 0.2,
00223 0.01,
00224 };
00225 static Double_t triggerParm_MC[]={8.0,
00226 4.0,
00227 0.95
00228 };
00229
00230
00231
00232 static Double_t triggerParm_MCSYS1[]={7.7,
00233 3.7,
00234 0.95
00235 };
00236
00237 static Double_t parm_MC_SYS2[]={
00238 0.0,
00239 0.0,
00240 3.5,
00241 100.,
00242 -2.0,
00243 0.1,
00244 -0.5,
00245 0.2,
00246 0.01,
00247 };
00248
00249
00250
00251
00252 TH1D *histDetMC[3];
00253 histDetMC[0]=new TH1D("PT(MC,signal)",";Pt(obs)",nDet,xminDet,xmaxDet);
00254 histDetMC[1]=new TH1D("PT(MC,bgr1)",";Pt(obs)",nDet,xminDet,xmaxDet);
00255 histDetMC[2]=new TH1D("PT(MC,bgr2)",";Pt(obs)",nDet,xminDet,xmaxDet);
00256 TH1D *histDetMCall=new TH1D("PT(MC,all)",";Pt(obs)",nDet,xminDet,xmaxDet);
00257 TH1D *histDetMCbgr=new TH1D("PT(MC,bgr)",";Pt(obs)",nDet,xminDet,xmaxDet);
00258
00259
00260 TH1D *histDetMCbbb[3];
00261 histDetMCbbb[0]=new TH1D("PT(MC,sig)",";Pt(obs,bbb)",nGen,xminGen,xmaxGen);
00262 histDetMCbbb[1]=new TH1D("PT(MC,bgr1)bbb",";Pt(obs,bbb)",nGen,xminGen,xmaxGen);
00263 histDetMCbbb[2]=new TH1D("PT(MC,bgr2)bbb",";Pt(obs,bbb)",nGen,xminGen,xmaxGen);
00264
00265 TH1D *histGenMC=new TH1D("PT(MC,signal,gen)bbb",";Pt(gen)",
00266 nGen,xminGen,xmaxGen);
00267 TH2D *histGenDet=new TH2D("PT(MC,signal,gen,det)",";Pt(gen);Pt(det)",
00268 nGen,xminGen,xmaxGen,nDet,xminDet,xmaxDet);
00269
00270
00271 Double_t lumiWeight = lumiData/lumiMC;
00272 intLumi=0.0;
00273 while(intLumi<lumiMC) {
00274 Int_t iTypeGen;
00275 Bool_t isTriggered;
00276 Double_t ptGen;
00277 Double_t ptObs=GenerateEvent(parm_MC,triggerParm_MC,&intLumi,&isTriggered,
00278 &ptGen,&iTypeGen);
00279
00280
00281 if(isTriggered) {
00282
00283 histDetMC[iTypeGen]->Fill(ptObs,lumiWeight);
00284 histDetMCbbb[iTypeGen]->Fill(ptObs,lumiWeight);
00285
00286 histDetMCall->Fill(ptObs,lumiWeight);
00287
00288 if(iTypeGen>0) {
00289 histDetMCbgr->Fill(ptObs,lumiWeight);
00290 }
00291 }
00292
00293 if(iTypeGen==0) {
00294 histGenMC->Fill(ptGen,lumiWeight);
00295 }
00296
00297 if(iTypeGen==0) {
00298
00299 if(isTriggered) {
00300 histGenDet->Fill(ptGen,ptObs,lumiWeight);
00301 } else {
00302
00303 histGenDet->Fill(ptGen,-999.,lumiWeight);
00304 }
00305 }
00306 }
00307
00308
00309
00310 TH2D *histGenDetSYS1=new TH2D("PT(MC,signal,gen,det,sys1)",
00311 ";Pt(gen);Pt(obs)",
00312 nGen,xminGen,xmaxGen,nDet,xminDet,xmaxDet);
00313 TH1D *histGenSYS1=new TH1D("PT(MC,signal,gen,sys1)",
00314 ";Pt(obs)",
00315 nGen,xminGen,xmaxGen);
00316 TH1D *histDetSYS1bbb=new TH1D("PT(MC,signal,det,sys1,bbb)",
00317 ";Pt(obs)",
00318 nGen,xminGen,xmaxGen);
00319 intLumi=0.;
00320 while(intLumi<lumiMC) {
00321 Int_t iTypeGen;
00322 Bool_t isTriggered;
00323 Double_t ptGen;
00324 Double_t ptObs=GenerateEvent(parm_MC,triggerParm_MCSYS1,&intLumi,
00325 &isTriggered,&ptGen,&iTypeGen);
00326
00327 if(iTypeGen==0) {
00328
00329 if(isTriggered) {
00330 histGenDetSYS1->Fill(ptGen,ptObs,lumiWeight);
00331 } else {
00332
00333 histGenDetSYS1->Fill(ptGen,-999.,lumiWeight);
00334 }
00335 }
00336
00337 if(iTypeGen==0) {
00338 histGenSYS1->Fill(ptGen,lumiWeight);
00339
00340 if(isTriggered) {
00341 histDetSYS1bbb->Fill(ptObs,lumiWeight);
00342 }
00343 }
00344 }
00345
00346
00347 TH2D *histGenDetSYS2=new TH2D("PT(MC,signal,gen,det,sys2)",
00348 ";Pt(gen);Pt(det)",
00349 nGen,xminGen,xmaxGen,nDet,xminDet,xmaxDet);
00350 TH1D *histGenSYS2=new TH1D("PT(MC,signal,gen,sys2)",
00351 ";Pt(obs)",
00352 nGen,xminGen,xmaxGen);
00353 TH1D *histDetSYS2bbb=new TH1D("PT(MC,signal,det,sys2,bbb)",
00354 ";Pt(obs)",
00355 nGen,xminGen,xmaxGen);
00356 intLumi=0.0;
00357 while(intLumi<lumiMC) {
00358 Int_t iTypeGen;
00359 Bool_t isTriggered;
00360 Double_t ptGen;
00361 Double_t ptObs=GenerateEvent(parm_MC_SYS2,triggerParm_MC,&intLumi,
00362 &isTriggered,&ptGen,&iTypeGen);
00363
00364 if(iTypeGen==0) {
00365
00366 if(isTriggered) {
00367 histGenDetSYS2->Fill(ptGen,ptObs,lumiWeight);
00368 } else {
00369
00370 histGenDetSYS2->Fill(ptGen,-999.,lumiWeight);
00371 }
00372 }
00373 if(iTypeGen==0) {
00374
00375 histGenSYS2->Fill(ptGen,lumiWeight);
00376
00377 if(isTriggered) {
00378 histDetSYS2bbb->Fill(ptObs,lumiWeight);
00379 }
00380 }
00381 }
00382
00383
00384
00385
00386
00387
00388
00389
00390
00391
00392
00393
00394
00395
00396
00397
00398
00399
00400
00401
00402
00403 TUnfoldSys unfold(histGenDet,TUnfold::kHistMapOutputHoriz,
00404 TUnfold::kRegModeSize,TUnfold::kEConstraintNone);
00405
00406
00407 unfold.SetInput(histDetDATA);
00408
00409
00410
00411 Double_t scale_bgr=1.0;
00412 Double_t dscale_bgr=0.2;
00413 unfold.SubtractBackground(histDetMC[1],"background1",scale_bgr,dscale_bgr);
00414 unfold.SubtractBackground(histDetMC[2],"background2",scale_bgr,dscale_bgr);
00415
00416
00417
00418 unfold.AddSysError(histGenDetSYS1,"trigger_SYS1",
00419 TUnfold::kHistMapOutputHoriz,
00420 TUnfoldSys::kSysErrModeMatrix);
00421
00422 unfold.AddSysError(histGenDetSYS2,"signalshape_SYS2",
00423 TUnfold::kHistMapOutputHoriz,
00424 TUnfoldSys::kSysErrModeMatrix);
00425
00426
00427 Int_t nScan=30;
00428 TSpline *logTauX,*logTauY;
00429 TGraph *lCurve;
00430
00431
00432 Int_t iBest=unfold.ScanLcurve(nScan,0.,0.,&lCurve,&logTauX,&logTauY);
00433
00434
00435 Double_t t[1],x[1],y[1];
00436 logTauX->GetKnot(iBest,t[0],x[0]);
00437 logTauY->GetKnot(iBest,t[0],y[0]);
00438 TGraph *bestLcurve=new TGraph(1,x,y);
00439 TGraph *bestLogTauLogChi2=new TGraph(1,t,x);
00440
00441
00442
00443
00444
00445
00446
00447
00448 Int_t *binMap=new Int_t[nGen+2];
00449 binMap[0] = -1;
00450 binMap[nGen+1]=-1;
00451 for(Int_t i=1;i<=nGen;i++) {
00452 binMap[i]=i;
00453 }
00454
00455
00456
00457
00458
00459 TH1D *histUnfoldStatBgr=new TH1D("PT(unfold,signal,stat+bgr)",";Pt(gen)",
00460 nGen,xminGen,xmaxGen);
00461 unfold.GetOutput(histUnfoldStatBgr,binMap);
00462
00463
00464 TH2D *histEmatStat=new TH2D("ErrMat(stat)",";Pt(gen);Pt(gen)",
00465 nGen,xminGen,xmaxGen,nGen,xminGen,xmaxGen);
00466 unfold.GetEmatrixInput(histEmatStat,binMap);
00467
00468
00469
00470 TH2D *histEmatTotal=new TH2D("ErrMat(total)",";Pt(gen);Pt(gen)",
00471 nGen,xminGen,xmaxGen,nGen,xminGen,xmaxGen);
00472 unfold.GetEmatrixTotal(histEmatTotal,binMap);
00473
00474
00475
00476 TH1D *histUnfoldStat=new TH1D("PT(unfold,signal,stat)",";Pt(gen)",
00477 nGen,xminGen,xmaxGen);
00478 TH1D *histUnfoldTotal=new TH1D("PT(unfold,signal,total)",";Pt(gen)",
00479 nGen,xminGen,xmaxGen);
00480 for(Int_t i=0;i<nGen+2;i++) {
00481 Double_t c=histUnfoldStatBgr->GetBinContent(i);
00482
00483
00484 histUnfoldStat->SetBinContent(i,c);
00485 histUnfoldStat->SetBinError
00486 (i,TMath::Sqrt(histEmatStat->GetBinContent(i,i)));
00487
00488
00489 histUnfoldTotal->SetBinContent(i,c);
00490 histUnfoldTotal->SetBinError
00491 (i,TMath::Sqrt(histEmatTotal->GetBinContent(i,i)));
00492 }
00493
00494
00495 TH2D *histCorr=new TH2D("Corr(total)",";Pt(gen);Pt(gen)",
00496 nGen,xminGen,xmaxGen,nGen,xminGen,xmaxGen);
00497 for(Int_t i=0;i<nGen+2;i++) {
00498 Double_t ei,ej;
00499 ei=TMath::Sqrt(histEmatTotal->GetBinContent(i,i));
00500 if(ei<=0.0) continue;
00501 for(Int_t j=0;j<nGen+2;j++) {
00502 ej=TMath::Sqrt(histEmatTotal->GetBinContent(j,j));
00503 if(ej<=0.0) continue;
00504 histCorr->SetBinContent(i,j,histEmatTotal->GetBinContent(i,j)/ei/ej);
00505 }
00506 }
00507
00508 delete [] binMap;
00509
00510
00511
00512
00513 TCanvas output;
00514 output.Divide(3,2);
00515 output.cd(1);
00516 histDetDATA->SetMinimum(0.0);
00517 histDetDATA->Draw("E");
00518 histDetMCall->SetMinimum(0.0);
00519 histDetMCall->SetLineColor(kBlue);
00520 histDetMCbgr->SetLineColor(kRed);
00521 histDetMC[1]->SetLineColor(kCyan);
00522 histDetMCall->Draw("SAME HIST");
00523 histDetMCbgr->Draw("SAME HIST");
00524 histDetMC[1]->Draw("SAME HIST");
00525
00526 output.cd(2);
00527 histUnfoldTotal->SetMinimum(0.0);
00528 histUnfoldTotal->SetMaximum(histUnfoldTotal->GetMaximum()*1.5);
00529
00530 histUnfoldTotal->Draw("E");
00531
00532 histUnfoldStatBgr->Draw("SAME E1");
00533
00534 histUnfoldStat->Draw("SAME E1");
00535 histGenDATA->Draw("SAME HIST");
00536 histGenMC->SetLineColor(kBlue);
00537 histGenMC->Draw("SAME HIST");
00538
00539 output.cd(3);
00540 histGenDet->SetLineColor(kBlue);
00541 histGenDet->Draw("BOX");
00542
00543
00544 output.cd(4);
00545 logTauX->Draw();
00546 bestLogTauLogChi2->SetMarkerColor(kRed);
00547 bestLogTauLogChi2->Draw("*");
00548
00549
00550 output.cd(5);
00551 lCurve->Draw("AL");
00552 bestLcurve->SetMarkerColor(kRed);
00553 bestLcurve->Draw("*");
00554
00555
00556 output.cd(6);
00557 histCorr->Draw("BOX");
00558
00559 output.SaveAs("testUnfold3.ps");
00560
00561
00562
00563
00564 for(Int_t i=1;i<=nGen;i++) {
00565
00566 Double_t data=histDetDATAbbb->GetBinContent(i);
00567 Double_t errData=histDetDATAbbb->GetBinError(i);
00568
00569
00570 Double_t data_bgr=data;
00571 Double_t errData_bgr=errData;
00572 for(Int_t j=1;j<=2;j++) {
00573 data_bgr -= scale_bgr*histDetMCbbb[j]->GetBinContent(i);
00574 errData_bgr = TMath::Sqrt(errData_bgr*errData_bgr+
00575 scale_bgr*histDetMCbbb[j]->GetBinError(i)*
00576 scale_bgr*histDetMCbbb[j]->GetBinError(i)+
00577 dscale_bgr*histDetMCbbb[j]->GetBinContent(i)*
00578 dscale_bgr*histDetMCbbb[j]->GetBinContent(i)
00579 );
00580 }
00581
00582 Double_t fCorr=(histGenMC->GetBinContent(i)/
00583 histDetMCbbb[0]->GetBinContent(i));
00584 Double_t data_bbb= data_bgr *fCorr;
00585
00586 Double_t errData_stat_bbb = errData*fCorr;
00587
00588 Double_t errData_statbgr_bbb = errData_bgr*fCorr;
00589
00590
00591 Double_t fCorr_1=(histGenSYS1->GetBinContent(i)/
00592 histDetSYS1bbb->GetBinContent(i));
00593 Double_t shift_sys1= data_bgr*fCorr_1 - data_bbb;
00594 Double_t fCorr_2=(histGenSYS2->GetBinContent(i)/
00595 histDetSYS2bbb->GetBinContent(i));
00596 Double_t shift_sys2= data_bgr*fCorr_2 - data_bbb;
00597
00598 cout<<data_bbb<<" "<<shift_sys1<<" "<<shift_sys2<<"\n";
00599
00600
00601 Double_t errData_total_bbb=
00602 TMath::Sqrt(errData_statbgr_bbb*errData_statbgr_bbb
00603 +shift_sys1*shift_sys1
00604 +shift_sys2*shift_sys2);
00605
00606
00607 Double_t data_unfold= histUnfoldStat->GetBinContent(i);
00608 Double_t errData_stat_unfold=histUnfoldStat->GetBinError(i);
00609 Double_t errData_statbgr_unfold=histUnfoldStatBgr->GetBinError(i);
00610 Double_t errData_total_unfold=histUnfoldTotal->GetBinError(i);
00611
00612
00613 std::cout<<"Bin "<<i<<": true "<<histGenDATA->GetBinContent(i)
00614 <<" unfold: "<<data_unfold
00615 <<" +/- "<<errData_stat_unfold<<" (stat)"
00616 <<" +/- "<<TMath::Sqrt(errData_statbgr_unfold*
00617 errData_statbgr_unfold-
00618 errData_stat_unfold*
00619 errData_stat_unfold)<<" (bgr)"
00620 <<" +/- "<<TMath::Sqrt(errData_total_unfold*
00621 errData_total_unfold-
00622 errData_statbgr_unfold*
00623 errData_statbgr_unfold)<<" (sys)"<<"\n";
00624 std::cout<<"Bin "<<i<<": true "<<histGenDATA->GetBinContent(i)
00625 <<" binbybin: "<<data_bbb
00626 <<" +/- "<<errData_stat_bbb<<" (stat)"
00627 <<" +/- "<<TMath::Sqrt(errData_statbgr_bbb*
00628 errData_statbgr_bbb-
00629 errData_stat_bbb*
00630 errData_stat_bbb)<<" (bgr)"
00631 <<" +/- "<<TMath::Sqrt(errData_total_bbb*
00632 errData_total_bbb-
00633 errData_statbgr_bbb*
00634 errData_statbgr_bbb)<<" (sys)"
00635 <<"\n";
00636 }
00637 }