00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027 #include <TError.h>
00028 #include <TMath.h>
00029 #include <TCanvas.h>
00030 #include <TRandom3.h>
00031 #include <TFitter.h>
00032 #include <TF1.h>
00033 #include <TStyle.h>
00034 #include <TVector.h>
00035 #include <TGraph.h>
00036
00037 #include <TUnfoldSys.h>
00038
00039
00040
00041 using namespace std;
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
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 TRandom *rnd=0;
00085
00086 TH2D *gHistInvEMatrix;
00087
00088 TVirtualFitter *gFitter=0;
00089
00090 void chisquare_corr(Int_t &npar, Double_t * , Double_t &f, Double_t *u, Int_t ) {
00091
00092
00093
00094
00095 Double_t x;
00096 TH1 *hfit = (TH1*)gFitter->GetObjectFit();
00097 TF1 *f1 = (TF1*)gFitter->GetUserFunc();
00098
00099 f1->InitArgs(&x,u);
00100 npar = f1->GetNpar();
00101 f = 0;
00102
00103 Int_t npfit = 0;
00104 Int_t nPoints=hfit->GetNbinsX();
00105 Double_t *df=new Double_t[nPoints];
00106 for (Int_t i=0;i<nPoints;i++) {
00107 x = hfit->GetBinCenter(i+1);
00108 TF1::RejectPoint(kFALSE);
00109 df[i] = f1->EvalPar(&x,u)-hfit->GetBinContent(i+1);
00110 if (TF1::RejectedPoint()) df[i]=0.0;
00111 else npfit++;
00112 }
00113 for (Int_t i=0;i<nPoints;i++) {
00114 for (Int_t j=0;j<nPoints;j++) {
00115 f += df[i]*df[j]*gHistInvEMatrix->GetBinContent(i+1,j+1);
00116 }
00117 }
00118 delete[] df;
00119 f1->SetNumberFitPoints(npfit);
00120 }
00121
00122 Double_t bw_func(Double_t *x,Double_t *par) {
00123 Double_t dm=x[0]-par[1];
00124 return par[0]/(dm*dm+par[2]*par[2]);
00125 }
00126
00127
00128
00129
00130
00131
00132 Double_t GenerateEvent(Double_t bgr,
00133 Double_t mass,
00134 Double_t gamma)
00135 {
00136 Double_t t;
00137 if(rnd->Rndm()>bgr) {
00138
00139
00140 do {
00141 do {
00142 t=rnd->Rndm();
00143 } while(t>=1.0);
00144 t=TMath::Tan((t-0.5)*TMath::Pi())*gamma+mass;
00145 } while(t<=0.0);
00146 return t;
00147 } else {
00148
00149
00150
00151 static Double_t const E0=2.4;
00152 static Double_t const N0=2.9;
00153 do {
00154 do {
00155 t=rnd->Rndm();
00156 } while(t>=1.0);
00157
00158
00159 t= -(TMath::Power(1.-t,1./(1.-N0))-1.0)*E0;
00160 } while(t>=0.0);
00161 return t;
00162 }
00163 }
00164
00165
00166
00167
00168
00169
00170 Double_t DetectorEvent(Double_t mTrue) {
00171
00172 static Double_t frac=0.1;
00173 static Double_t wideBias=0.03;
00174 static Double_t wideSigma=0.5;
00175 static Double_t smallBias=0.0;
00176 static Double_t smallSigma=0.1;
00177 if(rnd->Rndm()>frac) {
00178 return rnd->Gaus(mTrue+smallBias,smallSigma);
00179 } else {
00180 return rnd->Gaus(mTrue+wideBias,wideSigma);
00181 }
00182 }
00183
00184 int testUnfold1()
00185 {
00186
00187 TH1::SetDefaultSumw2();
00188
00189
00190 gStyle->SetOptFit(1111);
00191
00192
00193 rnd=new TRandom3();
00194
00195
00196 Double_t const luminosityData=100000;
00197 Double_t const luminosityMC=1000000;
00198 Double_t const crossSection=1.0;
00199
00200 Int_t const nDet=250;
00201 Int_t const nGen=100;
00202 Double_t const xminDet=0.0;
00203 Double_t const xmaxDet=10.0;
00204 Double_t const xminGen=0.0;
00205 Double_t const xmaxGen=10.0;
00206
00207
00208
00209
00210 TH1D *histMgenMC=new TH1D("MgenMC",";mass(gen)",nGen,xminGen,xmaxGen);
00211 TH1D *histMdetMC=new TH1D("MdetMC",";mass(det)",nDet,xminDet,xmaxDet);
00212 TH2D *histMdetGenMC=new TH2D("MdetgenMC",";mass(det);mass(gen)",
00213 nDet,xminDet,xmaxDet,nGen,xminGen,xmaxGen);
00214 Int_t neventMC=rnd->Poisson(luminosityMC*crossSection);
00215 for(Int_t i=0;i<neventMC;i++) {
00216 Double_t mGen=GenerateEvent(0.3,
00217 4.0,
00218 0.2);
00219 Double_t mDet=DetectorEvent(TMath::Abs(mGen));
00220
00221
00222
00223
00224
00225
00226
00227
00228 histMgenMC->Fill(mGen,luminosityData/luminosityMC);
00229
00230 histMdetMC->Fill(mDet,luminosityData/luminosityMC);
00231
00232
00233
00234
00235
00236
00237
00238
00239
00240
00241
00242
00243
00244
00245 histMdetGenMC->Fill(mDet,mGen,luminosityData/luminosityMC);
00246 }
00247
00248
00249
00250
00251
00252 TH2D *histMdetGenSysMC=new TH2D("MdetgenSysMC",";mass(det);mass(gen)",
00253 nDet,xminDet,xmaxDet,nGen,xminGen,xmaxGen);
00254 neventMC=rnd->Poisson(luminosityMC*crossSection);
00255 for(Int_t i=0;i<neventMC;i++) {
00256 Double_t mGen=GenerateEvent
00257 (0.5,
00258 3.6,
00259 0.15);
00260 Double_t mDet=DetectorEvent(TMath::Abs(mGen));
00261 histMdetGenSysMC->Fill(mDet,mGen,luminosityData/luminosityMC);
00262 }
00263
00264
00265
00266
00267 TH1D *histMgenData=new TH1D("MgenData",";mass(gen)",nGen,xminGen,xmaxGen);
00268 TH1D *histMdetData=new TH1D("MdetData",";mass(det)",nDet,xminDet,xmaxDet);
00269 Int_t neventData=rnd->Poisson(luminosityData*crossSection);
00270 for(Int_t i=0;i<neventData;i++) {
00271 Double_t mGen=GenerateEvent(0.4,
00272 3.8,
00273 0.15);
00274 Double_t mDet=DetectorEvent(TMath::Abs(mGen));
00275
00276
00277 histMgenData->Fill(mGen);
00278
00279
00280 histMdetData->Fill(mDet);
00281 }
00282
00283
00284
00285
00286 TUnfoldSys unfold(histMdetGenMC,TUnfold::kHistMapOutputVert);
00287
00288
00289
00290
00291
00292
00293
00294 if(unfold.SetInput(histMdetData)>=10000) {
00295 std::cout<<"Unfolding result may be wrong\n";
00296 }
00297
00298
00299
00300
00301
00302 Int_t nScan=30;
00303
00304 Double_t tauMin=0.0;
00305 Double_t tauMax=0.0;
00306 Int_t iBest;
00307 TSpline *logTauX,*logTauY;
00308 TGraph *lCurve;
00309
00310
00311 #ifdef VERBOSE_LCURVE_SCAN
00312 Int_t oldinfo=gErrorIgnoreLevel;
00313 gErrorIgnoreLevel=kInfo;
00314 #endif
00315
00316
00317 iBest=unfold.ScanLcurve(nScan,tauMin,tauMax,&lCurve,&logTauX,&logTauY);
00318
00319
00320 #ifdef VERBOSE_LCURVE_SCAN
00321 gErrorIgnoreLevel=oldinfo;
00322 #endif
00323
00324
00325
00326
00327
00328 Double_t SYS_ERROR1_MSTART=6;
00329 Double_t SYS_ERROR1_SIZE=0.1;
00330 TH2D *histMdetGenSys1=new TH2D("Mdetgensys1",";mass(det);mass(gen)",
00331 nDet,xminDet,xmaxDet,nGen,xminGen,xmaxGen);
00332 for(Int_t i=0;i<=nDet+1;i++) {
00333 if(histMdetData->GetBinCenter(i)>=SYS_ERROR1_MSTART) {
00334 for(Int_t j=0;j<=nGen+1;j++) {
00335 histMdetGenSys1->SetBinContent(i,j,SYS_ERROR1_SIZE);
00336 }
00337 }
00338 }
00339 unfold.AddSysError(histMdetGenSysMC,"SYSERROR_MC",TUnfold::kHistMapOutputVert,
00340 TUnfoldSys::kSysErrModeMatrix);
00341 unfold.AddSysError(histMdetGenSys1,"SYSERROR1",TUnfold::kHistMapOutputVert,
00342 TUnfoldSys::kSysErrModeRelative);
00343
00344
00345
00346
00347 std::cout<<"tau="<<unfold.GetTau()<<"\n";
00348 std::cout<<"chi**2="<<unfold.GetChi2A()<<"+"<<unfold.GetChi2L()
00349 <<" / "<<unfold.GetNdf()<<"\n";
00350 std::cout<<"chi**2(sys)="<<unfold.GetChi2Sys()<<"\n";
00351
00352
00353
00354
00355
00356 Double_t t[1],x[1],y[1];
00357 logTauX->GetKnot(iBest,t[0],x[0]);
00358 logTauY->GetKnot(iBest,t[0],y[0]);
00359 TGraph *bestLcurve=new TGraph(1,x,y);
00360 TGraph *bestLogTauLogChi2=new TGraph(1,t,x);
00361
00362
00363
00364
00365
00366
00367
00368
00369
00370
00371
00372
00373
00374 Int_t *binMap=new Int_t[nGen+2];
00375 for(Int_t i=1;i<=nGen;i++) binMap[i]=i;
00376 binMap[0]=-1;
00377 binMap[nGen+1]=-1;
00378
00379
00380 TH1D *histMunfold=new TH1D("Unfolded",";mass(gen)",nGen,xminGen,xmaxGen);
00381 unfold.GetOutput(histMunfold,binMap);
00382
00383
00384 TH1D *histMdetFold=unfold.GetFoldedOutput("FoldedBack",";mass(det)",
00385 xminDet,xmaxDet);
00386
00387
00388
00389 TH2D *histEmatData=new TH2D("EmatData",";mass(gen);mass(gen)",
00390 nGen,xminGen,xmaxGen,nGen,xminGen,xmaxGen);
00391 unfold.GetEmatrix(histEmatData,binMap);
00392
00393
00394
00395
00396 TH2D *histEmatTotal=new TH2D("EmatTotal",";mass(gen);mass(gen)",
00397 nGen,xminGen,xmaxGen,nGen,xminGen,xmaxGen);
00398 unfold.GetEmatrixTotal(histEmatTotal,binMap);
00399
00400
00401 TH1D *histTotalError=
00402 new TH1D("TotalError",";mass(gen)",nGen,xminGen,xmaxGen);
00403 for(Int_t bin=1;bin<=nGen;bin++) {
00404 histTotalError->SetBinContent(bin,histMunfold->GetBinContent(bin));
00405 histTotalError->SetBinError
00406 (bin,TMath::Sqrt(histEmatTotal->GetBinContent(bin,bin)));
00407 }
00408
00409
00410
00411
00412 gHistInvEMatrix=new TH2D("invEmat",";mass(gen);mass(gen)",
00413 nGen,xminGen,xmaxGen,nGen,xminGen,xmaxGen);
00414 TH1D *histRhoi=new TH1D("rho_I",";mass(gen)",nGen,xminGen,xmaxGen);
00415 unfold.GetRhoI(histRhoi,gHistInvEMatrix,binMap);
00416
00417
00418 delete[] binMap;
00419 binMap=0;
00420
00421
00422
00423
00424
00425
00426
00427
00428 gFitter=TVirtualFitter::Fitter(histMunfold);
00429 gFitter->SetFCN(chisquare_corr);
00430
00431 TF1 *bw=new TF1("bw",bw_func,xminGen,xmaxGen,3);
00432 bw->SetParameter(0,1000.);
00433 bw->SetParameter(1,3.8);
00434 bw->SetParameter(2,0.2);
00435
00436
00437 histMunfold->Fit(bw,"UE");
00438
00439
00440
00441 TCanvas output;
00442 output.Divide(3,2);
00443
00444
00445
00446
00447
00448
00449
00450 output.cd(1);
00451 histMdetGenMC->Draw("BOX");
00452
00453
00454
00455
00456
00457 output.cd(2);
00458 histTotalError->SetLineColor(kBlue);
00459 histTotalError->Draw("E");
00460 histMunfold->SetLineColor(kGreen);
00461 histMunfold->Draw("SAME E1");
00462 histMgenData->SetLineColor(kRed);
00463 histMgenData->Draw("SAME");
00464 histMgenMC->Draw("SAME HIST");
00465
00466
00467
00468
00469
00470 output.cd(3);
00471 histMdetFold->SetLineColor(kBlue);
00472 histMdetFold->Draw();
00473 histMdetMC->Draw("SAME HIST");
00474 TH1D *histInput=unfold.GetInput("Minput",";mass(det)",xminDet,xmaxDet);
00475 histInput->SetLineColor(kRed);
00476 histInput->Draw("SAME");
00477
00478
00479 output.cd(4);
00480 histRhoi->Draw();
00481
00482
00483 output.cd(5);
00484 logTauX->Draw();
00485 bestLogTauLogChi2->SetMarkerColor(kRed);
00486 bestLogTauLogChi2->Draw("*");
00487
00488
00489 output.cd(6);
00490 lCurve->Draw("AL");
00491 bestLcurve->SetMarkerColor(kRed);
00492 bestLcurve->Draw("*");
00493
00494 output.SaveAs("testUnfold1.ps");
00495
00496 return 0;
00497 }