63 Int_t GenerateGenEvent(Int_t nmax,
const Double_t *probability) {
68 Double_t f=rnd->Rndm();
70 while((r<nmax)&&(f>=probability[r])) {
77 Double_t GenerateRecEvent(
const Double_t *shapeParm) {
86 Double_t f=rnd->Rndm();
89 r=rnd->Gaus(shapeParm[1],shapeParm[2]);
91 r=rnd->Rndm()*(shapeParm[4]-shapeParm[3])+shapeParm[3];
96 void testUnfold4(
bool printInfo =
false)
100 if (!printInfo) gErrorIgnoreLevel = kWarning;
103 TH1::SetDefaultSumw2();
109 Double_t
const nData0= 500.0;
110 Double_t
const nMC0 = 50000.0;
115 Double_t
const xminDet=0.0;
116 Double_t
const xmaxDet=15.0;
120 Double_t
const xminGen=-0.5;
121 Double_t
const xmaxGen= 2.5;
125 static const Double_t genFrac[]={0.3,0.6,0.1};
128 static const Double_t genShape[][5]=
129 {{1.0,2.0,1.5,0.,15.},
130 {1.0,7.0,2.5,0.,15.},
131 {0.0,0.0,0.0,0.,15.}};
135 TH1D *histDetDATA=
new TH1D(
"Yrec",
";DATA(Yrec)",nDet,xminDet,xmaxDet);
139 TH2D *histGenDetMC=
new TH2D(
"Yrec%Xgen",
"MC(Xgen,Yrec)",
140 nGen,xminGen,xmaxGen,nDet,xminDet,xmaxDet);
142 TH1D *histUnfold=
new TH1D(
"Xgen",
";DATA(Xgen)",nGen,xminGen,xmaxGen);
144 TH1D **histPullNC=
new TH1D* [nGen];
145 TH1D **histPullArea=
new TH1D* [nGen];
146 for(
int i=0;i<nGen;i++) {
147 histPullNC[i]=
new TH1D(TString::Format(
"PullNC%d",i),
"pull",15,-3.,3.);
148 histPullArea[i]=
new TH1D(TString::Format(
"PullArea%d",i),
"pull",15,-3.,3.);
152 cout<<
"TUnfold version is "<<TUnfold::GetTUnfoldVersion()<<
"\n";
154 for(
int itoy=0;itoy<1000;itoy++) {
155 if(!(itoy %10)) cout<<
"toy iteration: "<<itoy<<
"\n";
156 histDetDATA->Reset();
157 histGenDetMC->Reset();
159 Int_t nData=rnd->Poisson(nData0);
160 for(Int_t i=0;i<nData;i++) {
161 Int_t iGen=GenerateGenEvent(nGen,genFrac);
162 Double_t yObs=GenerateRecEvent(genShape[iGen]);
163 histDetDATA->Fill(yObs);
166 Int_t nMC=rnd->Poisson(nMC0);
167 for(Int_t i=0;i<nMC;i++) {
168 Int_t iGen=GenerateGenEvent(nGen,genFrac);
169 Double_t yObs=GenerateRecEvent(genShape[iGen]);
170 histGenDetMC->Fill(iGen,yObs);
180 TUnfoldSys unfold(histGenDetMC,TUnfold::kHistMapOutputHoriz,
181 TUnfold::kRegModeSize,TUnfold::kEConstraintNone);
183 unfold.SetInput(histDetDATA,0.0,1.0);
186 unfold.ScanLcurve(50,0.,0.,0,0,0);
189 unfold.GetOutput(histUnfold);
191 for(
int i=0;i<nGen;i++) {
192 histPullNC[i]->Fill((histUnfold->GetBinContent(i+1)-genFrac[i]*nData0)/
193 histUnfold->GetBinError(i+1));
197 unfold.SetConstraint(TUnfold::kEConstraintArea);
200 unfold.ScanLcurve(50,0.,0.,0,0,0);
203 unfold.GetOutput(histUnfold);
205 for(
int i=0;i<nGen;i++) {
206 histPullArea[i]->Fill((histUnfold->GetBinContent(i+1)-genFrac[i]*nData0)/
207 histUnfold->GetBinError(i+1));
214 gStyle->SetOptFit(1111);
216 for(
int i=0;i<nGen;i++) {
218 histPullNC[i]->Fit(
"gaus");
219 histPullNC[i]->Draw();
221 for(
int i=0;i<nGen;i++) {
223 histPullArea[i]->Fit(
"gaus");
224 histPullArea[i]->Draw();
226 output.SaveAs(
"testUnfold4.ps");