%%cpp -d #include #include #include #include #include #include #include #include #include "TUnfold.h" TRandom *rnd=nullptr; %%cpp -d Double_t GenerateEvent(Double_t bgr, // relative fraction of background Double_t mass, // peak position Double_t gamma) // peak width { Double_t t; if(rnd->Rndm()>bgr) { // generate signal event // with positive mass do { do { t=rnd->Rndm(); } while(t>=1.0); t=TMath::Tan((t-0.5)*TMath::Pi())*gamma+mass; } while(t<=0.0); return t; } else { // generate background event // generate events following a power-law distribution // f(E) = K * TMath::power((E0+E),N0) static Double_t const E0=2.4; static Double_t const N0=2.9; do { do { t=rnd->Rndm(); } while(t>=1.0); // the mass is returned negative // In our example a convenient way to indicate it is a background event. t= -(TMath::Power(1.-t,1./(1.-N0))-1.0)*E0; } while(t>=0.0); return t; } } %%cpp -d Double_t DetectorEvent(Double_t mTrue) { // smear by double-gaussian static Double_t frac=0.1; static Double_t wideBias=0.03; static Double_t wideSigma=0.5; static Double_t smallBias=0.0; static Double_t smallSigma=0.1; if(rnd->Rndm()>frac) { return rnd->Gaus(mTrue+smallBias,smallSigma); } else { return rnd->Gaus(mTrue+wideBias,wideSigma); } } TH1::SetDefaultSumw2(); rnd=new TRandom3(); Double_t const luminosityData=100000; Double_t const luminosityMC=1000000; Double_t const crossSection=1.0; Int_t const nDet=250; Int_t const nGen=100; Double_t const xminDet=0.0; Double_t const xmaxDet=10.0; Double_t const xminGen=0.0; Double_t const xmaxGen=10.0; TH1D *histMgenMC=new TH1D("MgenMC",";mass(gen)",nGen,xminGen,xmaxGen); TH1D *histMdetMC=new TH1D("MdetMC",";mass(det)",nDet,xminDet,xmaxDet); TH2D *histMdetGenMC=new TH2D("MdetgenMC",";mass(det);mass(gen)",nDet,xminDet,xmaxDet, nGen,xminGen,xmaxGen); Int_t neventMC=rnd->Poisson(luminosityMC*crossSection); for(Int_t i=0;iFill(mGen,luminosityData/luminosityMC); // reconstructed MC distribution (for comparison only) histMdetMC->Fill(mDet,luminosityData/luminosityMC); // matrix describing how the generator input migrates to the // reconstructed level. Unfolding input. // NOTE on underflow/overflow bins: // (1) the detector level under/overflow bins are used for // normalisation ("efficiency" correction) // in our toy example, these bins are populated from tails // of the initial MC distribution. // (2) the generator level underflow/overflow bins are // unfolded. In this example: // underflow bin: background events reconstructed in the detector // overflow bin: signal events generated at masses > xmaxDet // for the unfolded result these bins will be filled // -> the background normalisation will be contained in the underflow bin histMdetGenMC->Fill(mDet,mGen,luminosityData/luminosityMC); } TH1D *histMgenData=new TH1D("MgenData",";mass(gen)",nGen,xminGen,xmaxGen); TH1D *histMdetData=new TH1D("MdetData",";mass(det)",nDet,xminDet,xmaxDet); Int_t neventData=rnd->Poisson(luminosityData*crossSection); for(Int_t i=0;iFill(mGen); // reconstructed mass, unfolding input histMdetData->Fill(mDet); } TUnfold unfold(histMdetGenMC,TUnfold::kHistMapOutputVert, TUnfold::kRegModeNone); Double_t estimatedPeakPosition=3.8; Int_t nPeek=3; TUnfold::ERegMode regMode=TUnfold::kRegModeCurvature; Int_t iPeek=(Int_t)(nGen*(estimatedPeakPosition-xminGen)/(xmaxGen-xminGen) // offset 1.5 // accounts for start bin 1 // and rounding errors +0.5 +1.5); unfold.RegularizeBins(1,1,iPeek-nPeek,regMode); unfold.RegularizeBins(iPeek+nPeek,1,nGen-(iPeek+nPeek),regMode); if(unfold.SetInput(histMdetData,0.0)>=10000) { std::cout<<"Unfolding result may be wrong\n"; } Double_t tauMin=0.0; Double_t tauMax=0.0; Int_t nScan=30; Int_t iBest; TSpline *logTauX,*logTauY; TGraph *lCurve; iBest=unfold.ScanLcurve(nScan,tauMin,tauMax,&lCurve,&logTauX,&logTauY); std::cout<<"tau="<GetKnot(iBest,t[0],x[0]); logTauY->GetKnot(iBest,t[0],y[0]); TGraph *bestLcurve=new TGraph(1,x,y); TGraph *bestLogTauX=new TGraph(1,t,x); Int_t *binMap=new Int_t[nGen+2]; for(Int_t i=1;i<=nGen;i++) binMap[i]=i; binMap[0]=-1; binMap[nGen+1]=-1; TH1D *histMunfold=new TH1D("Unfolded",";mass(gen)",nGen,xminGen,xmaxGen); unfold.GetOutput(histMunfold,binMap); TH1D *histMdetFold=new TH1D("FoldedBack","mass(det)",nDet,xminDet,xmaxDet); unfold.GetFoldedOutput(histMdetFold); TH1D *histRhoi=new TH1D("rho_I","mass",nGen,xminGen,xmaxGen); unfold.GetRhoI(histRhoi,binMap); delete[] binMap; binMap=nullptr; TCanvas output; output.Divide(3,2); output.cd(1); histMdetGenMC->Draw("BOX"); output.cd(2); histMunfold->SetLineColor(kBlue); histMunfold->Draw(); histMgenData->SetLineColor(kRed); histMgenData->Draw("SAME"); histMgenMC->Draw("SAME HIST"); output.cd(3); histMdetFold->SetLineColor(kBlue); histMdetFold->Draw(); histMdetData->SetLineColor(kRed); histMdetData->Draw("SAME"); histMdetMC->Draw("SAME HIST"); output.cd(4); histRhoi->Draw(); output.cd(5); logTauX->Draw(); bestLogTauX->SetMarkerColor(kRed); bestLogTauX->Draw("*"); output.cd(6); lCurve->Draw("AL"); bestLcurve->SetMarkerColor(kRed); bestLcurve->Draw("*"); output.SaveAs("testUnfold2.ps"); return 0; gROOT->GetListOfCanvases()->Draw()